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
/
response_matrix.fpp
1514 lines (1270 loc) · 52.4 KB
/
response_matrix.fpp
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 response_matrix
use netcdf
use mpi
implicit none
public :: init_response_matrix, finish_response_matrix
public :: read_response_matrix
public :: response_matrix_initialized
private
logical :: response_matrix_initialized = .false.
integer, parameter :: mat_unit = 70
#if defined MPI && defined ISO_C_BINDING
integer :: window = MPI_WIN_NULL
#endif
contains
subroutine init_response_matrix
use linear_solve, only: lu_decomposition
use fields_arrays, only: response_matrix
use stella_layouts, only: vmu_lo
use stella_layouts, only: iv_idx, is_idx
use kt_grids, only: naky
use extended_zgrid, only: iz_low, iz_up
use extended_zgrid, only: neigen, ikxmod
use extended_zgrid, only: nsegments
use extended_zgrid, only: nzed_segment
use extended_zgrid, only: periodic
use job_manage, only: time_message
use mp, only: proc0, job, mp_abort
use run_parameters, only: mat_gen, lu_option_switch
use run_parameters, only: lu_option_none, lu_option_local, lu_option_global
use system_fortran, only: systemf
#if defined MPI && defined ISO_C_BINDING
use, intrinsic :: iso_c_binding, only: c_ptr, c_f_pointer
use mp, only: curr_focus, sgproc0, mp_comm, sharedsubprocs, scope, barrier
use mpi
#endif
implicit none
integer :: iky, ie, iseg, iz
integer :: ikx
integer :: nz_ext, nresponse
integer :: idx
integer :: izl_offset, izup
#if defined MPI && defined ISO_C_BINDING
integer :: prior_focus, ierr, nbytes_real, rec_len
integer :: disp_unit = 1
integer (kind=MPI_ADDRESS_KIND) :: win_size, real_size
integer*8 :: cur_pos
type(c_ptr) :: bptr, cptr
#endif
real :: dum
complex, dimension (:), allocatable :: phiext
complex, dimension (:,:), allocatable :: gext
logical :: debug = .false.
character(100) :: message_dgdphi, message_QN, message_lu
real, dimension (2) :: time_response_matrix_dgdphi
real, dimension (2) :: time_response_matrix_QN
real, dimension (2) :: time_response_matrix_lu
! Related to the saving of the the matrices in netcdf format
character(len=15) :: fmt, job_str
character(len=100) :: file_name
integer :: istatus
istatus = 0
if (proc0.and.debug) then
write (*,*) " "
write (*,'(A)') " ############################################################"
write (*,'(A)') " RESPONSE MATRIX"
write (*,'(A)') " ############################################################"
end if
message_dgdphi = ' calculate dgdphi: '
message_QN = ' calculate QN: '
message_lu = ' calculate LU: '
time_response_matrix_dgdphi = 0
time_response_matrix_QN = 0
time_response_matrix_lu = 0
! All matrices handled by processor i_proc and job are stored
! on a single file named: response_mat_job.iproc
fmt = '(I5.5)'
if (proc0.and.mat_gen) then
call systemf('mkdir -p mat')
write (job_str, '(I1.1)') job
file_name = './mat/response_mat_'//trim(job_str)
open(unit=mat_unit, status='replace', file=file_name, &
position='rewind', action='write', form='unformatted')
write(unit=mat_unit) naky
end if
if (response_matrix_initialized) return
response_matrix_initialized = .true.
if (.not.allocated(response_matrix)) allocate (response_matrix(naky))
#if defined ISO_C_BINDING && defined MPI
inquire(iolength=rec_len) dum !this will return 2 or 8 for double
if(rec_len.eq.1.or.rec_len.eq.4) then
nbytes_real = 4
real_size = 4_MPI_ADDRESS_KIND
else if (rec_len.eq.2.or.rec_len.eq.8) then
nbytes_real = 8
real_size = 8_MPI_ADDRESS_KIND
else
call mp_abort('failure retrieving the size of a real')
endif
! Create a single shared memory window for all the response matrices and
! permutation arrays.
! Creating a window for each matrix/array would lead to performance
! degradation on some clusters
if(window.eq.MPI_WIN_NULL) then
prior_focus = curr_focus
call scope(sharedsubprocs)
win_size = 0
if(sgproc0) then
do iky = 1, naky
do ie = 1, neigen(iky)
if (periodic(iky)) then
nresponse = nsegments(ie,iky)*nzed_segment
else
nresponse = nsegments(ie,iky)*nzed_segment+1
end if
win_size = win_size &
+ int(nresponse,MPI_ADDRESS_KIND)*4_MPI_ADDRESS_KIND &
+ int(nresponse**2,MPI_ADDRESS_KIND)*2*real_size
enddo
enddo
endif
call mpi_win_allocate_shared(win_size,disp_unit,MPI_INFO_NULL,mp_comm, &
bptr,window,ierr)
if(.not.sgproc0) then
!make sure all the procs have the right memory address
call mpi_win_shared_query(window, 0, win_size, disp_unit, bptr, ierr)
endif
call mpi_win_fence(0,window,ierr)
!the following is a hack that allows us to perform pointer arithmetic in Fortran
cur_pos = transfer(bptr,cur_pos)
call scope(prior_focus)
endif
#endif
! for a given ky and set of connected kx values
! give a unit impulse to phi at each zed location
! in the extended domain and solve for h(zed_extended,(vpa,mu,s))
do iky = 1, naky
if (proc0.and.mat_gen) THEN
write(unit=mat_unit) iky, neigen(iky)
end if
! the response matrix for each ky has neigen(ky)
! independent sets of connected kx values
if (.not.associated(response_matrix(iky)%eigen)) &
allocate (response_matrix(iky)%eigen(neigen(iky)))
if(proc0.and.debug) then
call time_message(.false., time_response_matrix_dgdphi, message_dgdphi)
end if
! loop over the sets of connected kx values
do ie = 1, neigen(iky)
! number of zeds x number of segments
nz_ext = nsegments(ie,iky)*nzed_segment+1
! treat zonal mode specially to avoid double counting
! as it is periodic
if (periodic(iky)) then
nresponse = nz_ext-1
else
nresponse = nz_ext
end if
if (proc0.and.mat_gen) then
write(unit=mat_unit) ie, nresponse
end if
#if !defined ISO_C_BINDING || !defined MPI
! for each ky and set of connected kx values,
! must have a response matrix that is N x N
! with N = number of zeds per 2pi segment x number of 2pi segments
if (.not.associated(response_matrix(iky)%eigen(ie)%zloc)) &
allocate (response_matrix(iky)%eigen(ie)%zloc(nresponse,nresponse))
! response_matrix%idx is needed to keep track of permutations
! to the response matrix made during LU decomposition
! it will be input to LU back substitution during linear solve
if (.not.associated(response_matrix(iky)%eigen(ie)%idx)) &
allocate (response_matrix(iky)%eigen(ie)%idx(nresponse))
#else
!exploit MPIs shared memory framework to reduce memory consumption of
!response matrices
if (.not.associated(response_matrix(iky)%eigen(ie)%zloc)) then
cptr = transfer(cur_pos,cptr)
call c_f_pointer(cptr,response_matrix(iky)%eigen(ie)%zloc,(/nresponse,nresponse/))
cur_pos = cur_pos + nresponse**2*2*nbytes_real
endif
if (.not.associated(response_matrix(iky)%eigen(ie)%idx)) then
cptr = transfer(cur_pos,cptr)
call c_f_pointer(cptr,response_matrix(iky)%eigen(ie)%idx,(/nresponse/))
cur_pos = cur_pos + nresponse*4
endif
#endif
allocate (gext(nz_ext,vmu_lo%llim_proc:vmu_lo%ulim_alloc))
allocate (phiext(nz_ext))
! idx is the index in the extended zed domain
! that we are giving a unit impulse
idx = 0
! loop over segments, starting with 1
! first segment is special because it has
! one more unique zed value than all others
! since domain is [z0-pi:z0+pi], including both endpoints
! i.e., one endpoint is shared with the previous segment
iseg = 1
! ikxmod gives the kx corresponding to iseg,ie,iky
ikx = ikxmod(iseg,ie,iky)
izl_offset = 0
! avoid double-counting of periodic points for zonal mode (and other periodic modes)
if (periodic(iky)) then
izup = iz_up(iseg)-1
else
izup = iz_up(iseg)
end if
! no need to obtain response to impulses at negative kx values
do iz = iz_low(iseg), izup
idx = idx + 1
call get_dgdphi_matrix_column (iky, ikx, iz, ie, idx, nz_ext, nresponse, phiext, gext)
end do
! once we have used one segments, remaining segments
! have one fewer unique zed point
izl_offset = 1
if (nsegments(ie,iky) > 1) then
do iseg = 2, nsegments(ie,iky)
ikx = ikxmod(iseg,ie,iky)
do iz = iz_low(iseg)+izl_offset, iz_up(iseg)
idx = idx + 1
call get_dgdphi_matrix_column (iky, ikx, iz, ie, idx, nz_ext, nresponse, phiext, gext)
end do
if (izl_offset == 0) izl_offset = 1
end do
end if
deallocate (gext,phiext)
enddo
!DSO - This ends parallelization over velocity space.
! At this point every processor has int dv dgdphi for a given ky
! and so the quasineutrality solve and LU decomposition can be
! parallelized locally if need be.
! This is preferable to parallelization over ky as the LU
! decomposition (and perhaps QN) will be dominated by the
! ky with the most connections
if(proc0.and.debug) then
call time_message(.true. , time_response_matrix_dgdphi, message_dgdphi)
call time_message(.false., time_response_matrix_QN, message_QN)
end if
#ifdef ISO_C_BINDING
call mpi_win_fence(0,window,ierr)
#endif
! solve quasineutrality
! for local stella, this is a diagonal process, but global stella
! may require something more sophisticated
! loop over the sets of connected kx values
do ie = 1, neigen(iky)
#if defined ISO_C_BINDING && defined MPI
if(sgproc0) then
#endif
! number of zeds x number of segments
nz_ext = nsegments(ie,iky)*nzed_segment+1
! treat zonal mode specially to avoid double counting
! as it is periodic
if (periodic(iky)) then
nresponse = nz_ext-1
else
nresponse = nz_ext
end if
allocate (phiext(nz_ext))
do idx = 1, nresponse
phiext(nz_ext) = 0.0
phiext(:nresponse) = response_matrix(iky)%eigen(ie)%zloc(:,idx)
call get_fields_for_response_matrix (phiext, iky, ie)
! next need to create column in response matrix from phiext
! negative sign because matrix to be inverted in streaming equation
! is identity matrix - response matrix
! add in contribution from identity matrix
phiext(idx) = phiext(idx)-1.0
response_matrix(iky)%eigen(ie)%zloc(:,idx) = -phiext(:nresponse)
end do
deallocate (phiext)
#if defined ISO_C_BINDING && defined MPI
endif
#endif
enddo
#ifdef ISO_C_BINDING
call mpi_win_fence(0,window,ierr)
#endif
if(proc0.and.debug) then
call time_message(.true. , time_response_matrix_QN, message_QN)
call time_message(.false., time_response_matrix_lu, message_lu)
end if
!now we have the full response matrix. Finally, perform its LU decomposition
#ifdef MPI
select case (lu_option_switch)
case (lu_option_global)
call parallel_LU_decomposition_global(iky)
case (lu_option_local)
#ifdef ISO_C_BINDING
call parallel_LU_decomposition_local(iky)
#else
call mp_abort('Stella must be built with HAS_ISO_BINDING in order to use local parallel LU decomposition.')
#endif
case default
#endif
do ie = 1, neigen(iky)
#ifdef ISO_C_BINDING
if(sgproc0) then
#endif
! now that we have the reponse matrix for this ky and set of connected kx values
!get the LU decomposition so we are ready to solve the linear system
call lu_decomposition (response_matrix(iky)%eigen(ie)%zloc, &
response_matrix(iky)%eigen(ie)%idx,dum)
#ifdef ISO_C_BINDING
endif
#endif
enddo
#ifdef MPI
end select
#endif
if(proc0.and.debug) then
call time_message(.true., time_response_matrix_lu, message_lu)
end if
time_response_matrix_dgdphi = 0
time_response_matrix_QN = 0
time_response_matrix_lu = 0
do ie = 1, neigen(iky)
if (proc0.and.mat_gen) then
write(unit=mat_unit) response_matrix(iky)%eigen(ie)%idx
write(unit=mat_unit) response_matrix(iky)%eigen(ie)%zloc
end if
end do
!if(proc0) write (*,*) 'job', iky, iproc, response_matrix(iky)%eigen(1)%zloc(5,:)
end do
#ifdef ISO_C_BINDING
call mpi_win_fence(0,window,ierr)
#endif
if (proc0.and.mat_gen) then
close(unit=mat_unit)
end if
if (proc0.and.debug) then
write (*,'(A)') " ############################################################"
write (*,'(A)') " "
end if
end subroutine init_response_matrix
subroutine read_response_matrix
use fields_arrays, only: response_matrix
use common_types, only: response_matrix_type
use kt_grids, only: naky
use extended_zgrid, only: neigen
use extended_zgrid, only: nsegments
use extended_zgrid, only: nzed_segment
use extended_zgrid, only: periodic
use mp, only: proc0, job, broadcast, mp_abort
implicit none
integer :: iky, ie, nz_ext
integer :: iky_dump, neigen_dump, naky_dump, nresponse_dump
integer :: nresponse
character(len=15) :: fmt, job_str
character(len=100) :: file_name
integer :: istatus, ie_dump, istat
logical, parameter :: debug=.false.
istatus = 0
! All matrices handled for the job i_job are read
! from a single file named: responst_mat.ijob by that
! jobs root process
if(proc0) then
fmt = '(I5.5)'
write (job_str, '(I1.1)') job
file_name = './mat/response_mat.'//trim(job_str)
open(unit=mat_unit, status='old', file=file_name, &
action='read', form='unformatted', iostat=istat)
if (istat /= 0) then
print *, 'Error opening response_matrix by root processor for job ', job_str
end if
read(unit=mat_unit) naky_dump
if(naky.ne.naky_dump) call mp_abort('mismatch in naky and naky_dump')
endif
if (.not.allocated(response_matrix)) allocate (response_matrix(naky))
do iky = 1, naky
if(proc0) then
read(unit=mat_unit) iky_dump, neigen_dump
if(iky_dump.ne.iky.or.neigen_dump.ne.neigen(iky)) &
call mp_abort('mismatch in iky_dump/neigen_dump')
endif
if (.not.associated(response_matrix(iky)%eigen)) &
allocate (response_matrix(iky)%eigen(neigen(iky)))
! loop over the sets of connected kx values
do ie = 1, neigen(iky)
! number of zeds x number of segments
nz_ext = nsegments(ie,iky)*nzed_segment+1
! treat zonal mode specially to avoid double counting
! as it is periodic
if (periodic(iky)) then
nresponse = nz_ext-1
else
nresponse = nz_ext
end if
if(proc0) then
read(unit=mat_unit) ie_dump, nresponse_dump
if(ie_dump.ne.ie.or.nresponse.ne.nresponse_dump) &
call mp_abort('mismatch in ie/nresponse_dump')
endif
! for each ky and set of connected kx values,
! must have a response matrix that is N x N
! with N = number of zeds per 2pi segment x number of 2pi segments
if (.not.associated(response_matrix(iky)%eigen(ie)%zloc)) &
allocate (response_matrix(iky)%eigen(ie)%zloc(nresponse,nresponse))
! response_matrix%idx is needed to keep track of permutations
! to the response matrix made during LU decomposition
! it will be input to LU back substitution during linear solve
if (.not.associated(response_matrix(iky)%eigen(ie)%idx)) &
allocate (response_matrix(iky)%eigen(ie)%idx(nresponse))
if(proc0) then
read(unit=mat_unit) response_matrix(iky)%eigen(ie)%idx
read(unit=mat_unit) response_matrix(iky)%eigen(ie)%zloc
endif
call broadcast(response_matrix(iky)%eigen(ie)%idx)
call broadcast(response_matrix(iky)%eigen(ie)%zloc)
enddo
enddo
if (proc0) close (mat_unit)
if (debug) then
print *, 'File', file_name, ' successfully read by root proc for job: ', job_str
end if
end subroutine read_response_matrix
subroutine get_dgdphi_matrix_column (iky, ikx, iz, ie, idx, nz_ext, nresponse, phiext, gext)
use stella_layouts, only: vmu_lo
use stella_layouts, only: iv_idx, imu_idx, is_idx
use stella_time, only: code_dt
use zgrid, only: delzed, nzgrid
use extended_zgrid, only: periodic
use species, only: spec
use stella_geometry, only: gradpar, dbdzed
use vpamu_grids, only: vpa, mu
use vpamu_grids, only: maxwell_vpa, maxwell_mu, maxwell_fac
use fields_arrays, only: response_matrix
use gyro_averages, only: aj0x
use run_parameters, only: driftkinetic_implicit
use run_parameters, only: maxwellian_inside_zed_derivative
use parallel_streaming, only: stream_tridiagonal_solve
use parallel_streaming, only: stream_sign
use run_parameters, only: zed_upwind, time_upwind
#if defined ISO_C_BINDING && defined MPI
use mp, only: sgproc0
#endif
implicit none
integer, intent (in) :: iky, ikx, iz, ie, idx, nz_ext, nresponse
complex, dimension (:), intent (in out) :: phiext
complex, dimension (:,vmu_lo%llim_proc:), intent (in out) :: gext
integer :: ivmu, iv, imu, is, ia
integer :: izp, izm
real :: mu_dbdzed_p, mu_dbdzed_m
real :: fac, fac0, fac1, gyro_fac
ia = 1
if (.not.maxwellian_inside_zed_derivative) then
! get -vpa*b.gradz*Ze/T*F0*d<phi>/dz corresponding to unit impulse in phi
do ivmu = vmu_lo%llim_proc, vmu_lo%ulim_proc
! initialize g to zero everywhere along extended zed domain
gext(:,ivmu) = 0.0
iv = iv_idx(vmu_lo,ivmu)
imu = imu_idx(vmu_lo,ivmu)
is = is_idx(vmu_lo,ivmu)
! give unit impulse to phi at this zed location
! and compute -vpa*b.gradz*Ze/T*d<phi>/dz*F0 (RHS of streaming part of GKE)
! NB: assuming equal spacing in zed below
! here, fac = -dt*(1+alph_t)/2*vpa*Ze/T*F0*J0/dz
! b.gradz left out because needs to be centred in zed
if (driftkinetic_implicit) then
gyro_fac = 1.0
else
gyro_fac = aj0x(iky,ikx,iz,ivmu)
end if
! 0.125 to account for two linear interpolations
fac = -0.125*(1.+time_upwind)*code_dt*vpa(iv)*spec(is)%stm_psi0 &
*gyro_fac*spec(is)%zt/delzed(0)*maxwell_vpa(iv,is)*maxwell_fac(is)
! In the following, gradpar and maxwell_mu are interpolated separately
! to ensure consistency to what is done in parallel_streaming.f90
! stream_sign < 0 corresponds to positive advection speed
if (stream_sign(iv)<0) then
if (iz > -nzgrid) then
! fac0 is the factor multiplying delphi on the RHS
! of the homogeneous GKE at this zed index
fac0 = fac*((1.+zed_upwind)*gradpar(iz) &
+ (1.-zed_upwind)*gradpar(iz-1)) &
*((1.+zed_upwind)*maxwell_mu(ia,iz,imu,is) &
+ (1.-zed_upwind)*maxwell_mu(ia,iz-1,imu,is))
! fac1 is the factor multiplying delphi on the RHS
! of the homogeneous GKE at the zed index to the right of
! this one
if (iz < nzgrid) then
fac1 = fac*((1.+zed_upwind)*gradpar(iz+1) &
+ (1.-zed_upwind)*gradpar(iz)) &
*((1.+zed_upwind)*maxwell_mu(ia,iz+1,imu,is) &
+ (1.-zed_upwind)*maxwell_mu(ia,iz,imu,is))
else
fac1 = fac*((1.+zed_upwind)*gradpar(-nzgrid+1) &
+ (1.-zed_upwind)*gradpar(nzgrid)) &
*((1.+zed_upwind)*maxwell_mu(ia,-nzgrid+1,imu,is) &
+ (1.-zed_upwind)*maxwell_mu(ia,nzgrid,imu,is))
end if
else
! fac0 is the factor multiplying delphi on the RHS
! of the homogeneous GKE at this zed index
fac0 = fac*((1.+zed_upwind)*gradpar(iz) &
+ (1.-zed_upwind)*gradpar(nzgrid-1)) &
*((1.+zed_upwind)*maxwell_mu(ia,iz,imu,is) &
+ (1.-zed_upwind)*maxwell_mu(ia,nzgrid-1,imu,is))
! fac1 is the factor multiplying delphi on the RHS
! of the homogeneous GKE at the zed index to the right of
! this one
fac1 = fac*((1.+zed_upwind)*gradpar(iz+1) &
+ (1.-zed_upwind)*gradpar(iz)) &
*((1.+zed_upwind)*maxwell_mu(ia,iz+1,imu,is) &
+ (1.-zed_upwind)*maxwell_mu(ia,iz,imu,is))
end if
gext(idx,ivmu) = fac0
if (idx < nz_ext) gext(idx+1,ivmu) = -fac1
! zonal mode BC is periodic instead of zero, so must
! treat specially
if (periodic(iky)) then
if (idx == 1) then
gext(nz_ext,ivmu) = fac0
else if (idx == nz_ext-1) then
gext(1,ivmu) = -fac1
end if
end if
else
if (iz < nzgrid) then
! fac0 is the factor multiplying delphi on the RHS
! of the homogeneous GKE at this zed index
fac0 = fac*((1.+zed_upwind)*gradpar(iz) &
+ (1.-zed_upwind)*gradpar(iz+1)) &
*((1.+zed_upwind)*maxwell_mu(ia,iz,imu,is) &
+ (1.-zed_upwind)*maxwell_mu(ia,iz+1,imu,is))
! fac1 is the factor multiplying delphi on the RHS
! of the homogeneous GKE at the zed index to the left of
! this one
if (iz > -nzgrid) then
fac1 = fac*((1.+zed_upwind)*gradpar(iz-1) &
+ (1.-zed_upwind)*gradpar(iz)) &
*((1.+zed_upwind)*maxwell_mu(ia,iz-1,imu,is) &
+ (1.-zed_upwind)*maxwell_mu(ia,iz,imu,is))
else
fac1 = fac*((1.+zed_upwind)*gradpar(nzgrid-1) &
+ (1.-zed_upwind)*gradpar(iz)) &
*((1.+zed_upwind)*maxwell_mu(ia,nzgrid-1,imu,is) &
+ (1.-zed_upwind)*maxwell_mu(ia,iz,imu,is))
end if
else
! fac0 is the factor multiplying delphi on the RHS
! of the homogeneous GKE at this zed index
fac0 = fac*((1.+zed_upwind)*gradpar(iz) &
+ (1.-zed_upwind)*gradpar(-nzgrid+1)) &
*((1.+zed_upwind)*maxwell_mu(ia,iz,imu,is) &
+ (1.-zed_upwind)*maxwell_mu(ia,-nzgrid+1,imu,is))
! fac1 is the factor multiplying delphi on the RHS
! of the homogeneous GKE at the zed index to the left of
! this one
fac1 = fac*((1.+zed_upwind)*gradpar(iz-1) &
+ (1.-zed_upwind)*gradpar(iz)) &
*((1.+zed_upwind)*maxwell_mu(ia,iz-1,imu,is) &
+ (1.-zed_upwind)*maxwell_mu(ia,iz,imu,is))
end if
gext(idx,ivmu) = -fac0
if (idx > 1) gext(idx-1,ivmu) = fac1
! zonal mode BC is periodic instead of zero, so must
! treat specially
if (periodic(iky)) then
if (idx == 1) then
gext(nz_ext,ivmu) = -fac0
gext(nz_ext-1,ivmu) = fac1
else if (idx == 2) then
gext(nz_ext,ivmu) = fac1
end if
end if
end if
! hack for now (duplicates much of the effort from sweep_zed_zonal)
if (periodic(iky)) then
call sweep_zed_zonal_response (iv, is, stream_sign(iv), gext(:,ivmu))
else
! invert parallel streaming equation to get g^{n+1} on extended zed grid
! (I + (1+alph)/2*dt*vpa)*g_{inh}^{n+1} = RHS = gext
call stream_tridiagonal_solve (iky, ie, iv, is, gext(:,ivmu))
end if
end do
else
! get -vpa*b.gradz*Ze/T*F0*d<phi>/dz corresponding to unit impulse in phi
do ivmu = vmu_lo%llim_proc, vmu_lo%ulim_proc
! initialize g to zero everywhere along extended zed domain
gext(:,ivmu) = 0.0
iv = iv_idx(vmu_lo,ivmu)
imu = imu_idx(vmu_lo,ivmu)
is = is_idx(vmu_lo,ivmu)
! give unit impulse to phi at this zed location
! and compute -vpa*b.gradz*Ze/T*d<phi>/dz*F0 (RHS of streaming part of GKE)
! NB: assuming equal spacing in zed below
! here, fac = -dt*(1+alph_t)/2*vpa*Ze/T*F0*J0/dz
! b.gradz left out because needs to be centred in zed
if (driftkinetic_implicit) then
gyro_fac = 1.0
else
gyro_fac = aj0x(iky,ikx,iz,ivmu)
end if
fac = -0.25*(1.+time_upwind)*code_dt*vpa(iv)*spec(is)%stm_psi0 &
*gyro_fac*spec(is)%zt*maxwell_vpa(iv,is)*maxwell_mu(ia,iz,imu,is)*maxwell_fac(is)
mu_dbdzed_p = 1./delzed(0)+mu(imu)*dbdzed(ia,iz)*(1.+zed_upwind)
mu_dbdzed_m = 1./delzed(0)+mu(imu)*dbdzed(ia,iz)*(1.-zed_upwind)
! stream_sign < 0 corresponds to positive advection speed
if (stream_sign(iv)<0) then
if (iz > -nzgrid) then
! fac0 is the factor multiplying delphi on the RHS
! of the homogeneous GKE at this zed index
fac0 = fac*((1.+zed_upwind)*gradpar(iz) &
+ (1.-zed_upwind)*gradpar(iz-1))*mu_dbdzed_p
! fac1 is the factor multiplying delphi on the RHS
! of the homogeneous GKE at the zed index to the right of
! this one
if (iz < nzgrid) then
izp = iz+1
else
izp = -nzgrid+1
end if
fac1 = fac*((1.+zed_upwind)*gradpar(izp) &
+ (1.-zed_upwind)*gradpar(iz))*mu_dbdzed_m
else
! fac0 is the factor multiplying delphi on the RHS
! of the homogeneous GKE at this zed index
fac0 = fac*((1.+zed_upwind)*gradpar(iz) &
+ (1.-zed_upwind)*gradpar(nzgrid-1))*mu_dbdzed_p
! fac1 is the factor multiplying delphi on the RHS
! of the homogeneous GKE at the zed index to the right of
! this one
fac1 = fac*((1.+zed_upwind)*gradpar(iz+1) &
+ (1.-zed_upwind)*gradpar(iz))*mu_dbdzed_m
end if
gext(idx,ivmu) = fac0
if (idx < nz_ext) gext(idx+1,ivmu) = -fac1
! zonal mode BC is periodic instead of zero, so must
! treat specially
if (periodic(iky)) then
if (idx == 1) then
gext(nz_ext,ivmu) = fac0
else if (idx == nz_ext-1) then
gext(1,ivmu) = -fac1
end if
end if
else
if (iz < nzgrid) then
! fac0 is the factor multiplying delphi on the RHS
! of the homogeneous GKE at this zed index
fac0 = fac*((1.+zed_upwind)*gradpar(iz) &
+ (1.-zed_upwind)*gradpar(iz+1))*mu_dbdzed_p
! fac1 is the factor multiplying delphi on the RHS
! of the homogeneous GKE at the zed index to the left of
! this one
if (iz > -nzgrid) then
izm = iz-1
else
izm = nzgrid-1
end if
fac1 = fac*((1.+zed_upwind)*gradpar(izm) &
+ (1.-zed_upwind)*gradpar(iz))*mu_dbdzed_m
else
! fac0 is the factor multiplying delphi on the RHS
! of the homogeneous GKE at this zed index
fac0 = fac*((1.+zed_upwind)*gradpar(iz) &
+ (1.-zed_upwind)*gradpar(-nzgrid+1))*mu_dbdzed_p
! fac1 is the factor multiplying delphi on the RHS
! of the homogeneous GKE at the zed index to the left of
! this one
fac1 = fac*((1.+zed_upwind)*gradpar(iz-1) &
+ (1.-zed_upwind)*gradpar(iz))*mu_dbdzed_m
end if
gext(idx,ivmu) = -fac0
if (idx > 1) gext(idx-1,ivmu) = fac1
! zonal mode BC is periodic instead of zero, so must
! treat specially
if (periodic(iky)) then
if (idx == 1) then
gext(nz_ext,ivmu) = -fac0
gext(nz_ext-1,ivmu) = fac1
else if (idx == 2) then
gext(nz_ext,ivmu) = fac1
end if
end if
end if
! hack for now (duplicates much of the effort from sweep_zed_zonal)
if (periodic(iky)) then
call sweep_zed_zonal_response (iv, is, stream_sign(iv), gext(:,ivmu))
else
! invert parallel streaming equation to get g^{n+1} on extended zed grid
! (I + (1+alph)/2*dt*vpa)*g_{inh}^{n+1} = RHS = gext
call stream_tridiagonal_solve (iky, ie, iv, is, gext(:,ivmu))
end if
end do
end if
! we now have g on the extended zed domain at this ky and set of connected kx values
! corresponding to a unit impulse in phi at this location
! now integrate over velocities to get a square response matrix
! (this ends the parallelization over velocity space, so every core should have a
! copy of phiext)
call integrate_over_velocity (gext, phiext, iky, ie)
#if !defined ISO_C_BINDING || !defined MPI
response_matrix(iky)%eigen(ie)%zloc(:,idx) = phiext(:nresponse)
#else
if(sgproc0) response_matrix(iky)%eigen(ie)%zloc(:,idx) = phiext(:nresponse)
#endif
end subroutine get_dgdphi_matrix_column
! subroutine get_phi_matrix
! end subroutine get_phi_matrix
subroutine sweep_zed_zonal_response (iv, is, sgn, g)
use zgrid, only: nzgrid, delzed, nztot
use run_parameters, only: zed_upwind, time_upwind
use parallel_streaming, only: stream_c
implicit none
integer, intent (in) :: iv, is, sgn
complex, dimension (:), intent (in out) :: g
integer :: iz, iz1, iz2
real :: fac1, fac2
complex, dimension (:), allocatable :: gcf, gpi
allocate (gpi(-nzgrid:nzgrid))
allocate (gcf(-nzgrid:nzgrid))
! ky=0 is 2pi periodic (no extended zgrid)
! decompose into complementary function + particular integral
! zero BC for particular integral
! unit BC for complementary function (no source)
if (sgn < 0) then
iz1 = -nzgrid ; iz2 = nzgrid
else
iz1 = nzgrid ; iz2 = -nzgrid
end if
gpi(iz1) = 0. ; gcf(iz1) = 1.
do iz = iz1-sgn, iz2, -sgn
fac1 = 1.0+zed_upwind+sgn*(1.0+time_upwind)*stream_c(iz,iv,is)/delzed(0)
fac2 = 1.0-zed_upwind-sgn*(1.0+time_upwind)*stream_c(iz,iv,is)/delzed(0)
gpi(iz) = (-gpi(iz+sgn)*fac2 + 2.0*g(iz+nzgrid+1))/fac1
gcf(iz) = -gcf(iz+sgn)*fac2/fac1
end do
! g = g_PI + (g_PI(pi)/(1-g_CF(pi))) * g_CF
g = gpi + (spread(gpi(iz2),1,nztot)/(1.-gcf(iz2)))*gcf
deallocate (gpi, gcf)
end subroutine sweep_zed_zonal_response
subroutine integrate_over_velocity(g,phi,iky,ie)
use stella_layouts, only: vmu_lo
use species, only: nspec, spec
use extended_zgrid, only: iz_low, iz_up
use extended_zgrid, only: ikxmod
use extended_zgrid, only: nsegments
use kt_grids, only: zonal_mode, akx
use vpamu_grids, only: integrate_species
use gyro_averages, only: gyro_average
use mp, only: sum_allreduce
implicit none
complex, dimension (:,vmu_lo%llim_proc:), intent (in) :: g
complex, dimension (:), intent (out) :: phi
integer, intent (in) :: iky, ie
integer :: idx, iseg, ikx, iz, ia
integer :: izl_offset
real, dimension (nspec) :: wgt
complex, dimension (:), allocatable :: g0
ia = 1
allocate (g0(vmu_lo%llim_proc:vmu_lo%ulim_alloc))
wgt = spec%z*spec%dens_psi0
phi = 0.
idx = 0 ; izl_offset = 0
iseg = 1
ikx = ikxmod(iseg,ie,iky)
if(zonal_mode(iky).and.abs(akx(ikx)) < epsilon(0.)) then
phi(:) = 0.0
return
endif
do iz = iz_low(iseg), iz_up(iseg)
idx = idx + 1
call gyro_average (g(idx,:), iky, ikx, iz, g0)
call integrate_species (g0, iz, wgt, phi(idx),reduce_in=.false.)
end do
izl_offset = 1
if (nsegments(ie,iky) > 1) then
do iseg = 2, nsegments(ie,iky)
ikx = ikxmod(iseg,ie,iky)
do iz = iz_low(iseg)+izl_offset, iz_up(iseg)
idx = idx + 1
call gyro_average (g(idx,:), iky, ikx, iz, g0)
call integrate_species (g0, iz, wgt, phi(idx),reduce_in=.false.)
end do
if (izl_offset == 0) izl_offset = 1
end do
end if
call sum_allreduce(phi)
end subroutine integrate_over_velocity
subroutine get_fields_for_response_matrix (phi, iky, ie)
use stella_layouts, only: vmu_lo
use species, only: spec
use species, only: has_electron_species
use stella_geometry, only: dl_over_b
use extended_zgrid, only: iz_low, iz_up
use extended_zgrid, only: ikxmod
use extended_zgrid, only: nsegments
use kt_grids, only: zonal_mode, akx
use dist_fn, only: adiabatic_option_switch
use dist_fn, only: adiabatic_option_fieldlineavg
use fields, only: gamtot3
use fields_arrays, only: gamtot
implicit none
complex, dimension (:), intent (inout) :: phi
integer, intent (in) :: iky, ie
integer :: idx, iseg, ikx, iz, ia
integer :: izl_offset
complex, dimension (:), allocatable :: g0
complex :: tmp
ia = 1
allocate (g0(vmu_lo%llim_proc:vmu_lo%ulim_alloc))
idx = 0 ; izl_offset = 0
iseg = 1
ikx = ikxmod(iseg,ie,iky)
if(zonal_mode(iky).and.abs(akx(ikx)) < epsilon(0.)) then
phi(:) = 0.0
return
endif
do iz = iz_low(iseg), iz_up(iseg)
idx = idx + 1
phi(idx) = phi(idx)/gamtot(iky,ikx,iz)
end do
izl_offset = 1
if (nsegments(ie,iky) > 1) then
do iseg = 2, nsegments(ie,iky)
ikx = ikxmod(iseg,ie,iky)
do iz = iz_low(iseg)+izl_offset, iz_up(iseg)
idx = idx + 1
phi(idx) = phi(idx)/gamtot(iky,ikx,iz)
end do
if (izl_offset == 0) izl_offset = 1
end do
end if
if (.not.has_electron_species(spec) .and. &
adiabatic_option_switch == adiabatic_option_fieldlineavg) then
if (zonal_mode(iky)) then
! no connections for ky = 0
iseg = 1
tmp = sum(dl_over_b(ia,:)*phi)
phi = phi + tmp*gamtot3(ikxmod(1,ie,iky),:)
end if
end if
deallocate (g0)
end subroutine get_fields_for_response_matrix
subroutine finish_response_matrix
use fields_arrays, only: response_matrix
#if !defined ISO_C_BINDING
implicit none
#else
use mpi
implicit none
integer :: ierr
if (window.ne.MPI_WIN_NULL) call mpi_win_free(window,ierr)
#endif
if (allocated(response_matrix)) deallocate (response_matrix)
response_matrix_initialized = .false.
end subroutine finish_response_matrix
#ifdef MPI
!----------------------------------------------------------!
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! !
! !!!!!!!!!!!!!!PARALLEL LU DECOMPOSITIONS!!!!!!!!!!!!!!!! !