This repository has been archived by the owner on Jul 20, 2021. It is now read-only.
-
Notifications
You must be signed in to change notification settings - Fork 5
/
Copy pathdissipation.f90
2731 lines (2307 loc) · 130 KB
/
dissipation.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 dissipation
implicit none
public :: init_dissipation, finish_dissipation
public :: init_collisions, collisions_initialized
public :: include_collisions
public :: include_krook_operator, update_delay_krook
public :: remove_zero_projection, project_out_zero
public :: advance_collisions_explicit, advance_collisions_implicit
public :: time_collisions
public :: hyper_dissipation
public :: advance_hyper_dissipation
public :: add_krook_operator
public :: collisions_implicit
public :: delay_krook, int_krook, int_proj
public :: vpa_operator, mu_operator
public :: cfl_dt_vpadiff, cfl_dt_mudiff
private
logical :: include_collisions, vpa_operator, mu_operator
logical :: collisions_implicit, include_krook_operator
logical :: momentum_conservation, energy_conservation
logical :: hyper_dissipation, remove_zero_projection
logical :: use_physical_ksqr
logical :: krook_odd,exclude_boundary_regions
real :: D_hyper, nu_krook, delay_krook, int_krook, int_proj
real :: cfl_dt_vpadiff, cfl_dt_mudiff
integer:: ikxmax_source
character(30) :: collision_model
integer :: nresponse_vpa = 1
integer :: nresponse_mu = 1
real :: cfac
real, dimension (:,:), allocatable :: aa_vpa, bb_vpa, cc_vpa
real, dimension (:,:,:), allocatable :: aa_mu, cc_mu
real, dimension (:,:), allocatable :: bb_mu
complex, dimension (:,:,:), allocatable :: vpadiff_response
integer, dimension (:,:), allocatable :: vpadiff_idx
complex, dimension (:,:,:), allocatable :: mudiff_response
integer, dimension (:,:), allocatable :: mudiff_idx
complex, dimension (:,:,:,:,:), allocatable :: aa_blcs, cc_blcs
complex, dimension (:,:,:,:), allocatable :: bb_blcs
complex, dimension (:,:,:,:), allocatable :: dd_vecs
complex, dimension (:,:,:), allocatable :: cdiffmat_band
complex, dimension (:,:,:), allocatable :: blockmatrix
integer, dimension (:,:), allocatable :: ipiv
real, dimension (:,:,:), allocatable :: nus, nuD, nupa, nux, mw
integer :: info
logical :: collisions_initialized = .false.
real, dimension (2,2) :: time_collisions = 0.
contains
subroutine init_dissipation
use mp, only: proc0
use kt_grids, only: nakx
use zgrid, only: nzgrid, ntubes
use stella_layouts, only: vmu_lo
use dist_fn_arrays, only: g_krook, g_proj
implicit none
call read_parameters
if (include_collisions) then
call init_collisions
else
if (proc0) then
write (*,'(A)') "############################################################"
write (*,'(A)') " COLLISIONS"
write (*,'(A)') "############################################################"
write (*,*) 'Coll. model: None'
write (*,*)
end if
end if
if(include_krook_operator.and..not.allocated(g_krook)) then
allocate (g_krook(nakx,-nzgrid:nzgrid,ntubes,vmu_lo%llim_proc:vmu_lo%ulim_alloc))
g_krook = 0.
endif
if(remove_zero_projection.and..not.allocated(g_proj)) then
allocate (g_proj(nakx,-nzgrid:nzgrid,ntubes,vmu_lo%llim_proc:vmu_lo%ulim_alloc))
g_proj = 0.
endif
int_krook = 0.
int_proj = 0.
end subroutine init_dissipation
subroutine read_parameters
use file_utils, only: input_unit_exist
use physics_flags, only: full_flux_surface, radial_variation
use mp, only: proc0, broadcast
use kt_grids, only: ikx_max, periodic_variation
implicit none
namelist /dissipation/ collision_model, hyper_dissipation, D_hyper, &
include_collisions, collisions_implicit, &
momentum_conservation, energy_conservation, &
vpa_operator, mu_operator, include_krook_operator, &
nu_krook, delay_krook, remove_zero_projection, &
ikxmax_source, cfac, krook_odd, use_physical_ksqr, &
exclude_boundary_regions
integer :: in_file
logical :: dexist
if (proc0) then
include_collisions = .false.
include_krook_operator = .false.
collisions_implicit = .true.
collision_model = "dougherty"
momentum_conservation = .true.
energy_conservation = .true.
vpa_operator = .true.
mu_operator = .true.
hyper_dissipation = .false.
use_physical_ksqr = .not.(full_flux_surface.or.radial_variation)
exclude_boundary_regions = radial_variation.and..not.periodic_variation
remove_zero_projection = .false.
D_hyper = 0.05
nu_krook = 0.05
delay_krook =0.02
ikxmax_source = 1 ! kx=0
if(periodic_variation) ikxmax_source = 2 ! kx=0 and kx=1
krook_odd = .true. ! damp only the odd mode that can affect profiles
cfac = 1
in_file = input_unit_exist("dissipation", dexist)
if (dexist) read (unit=in_file, nml=dissipation)
end if
ikxmax_source = min(ikxmax_source,ikx_max)
call broadcast (include_collisions)
call broadcast (include_krook_operator)
call broadcast (collisions_implicit)
call broadcast (collision_model)
call broadcast (momentum_conservation)
call broadcast (energy_conservation)
call broadcast (vpa_operator)
call broadcast (mu_operator)
call broadcast (hyper_dissipation)
call broadcast (exclude_boundary_regions)
call broadcast (use_physical_ksqr)
call broadcast (D_hyper)
call broadcast (nu_krook)
call broadcast (delay_krook)
call broadcast (ikxmax_source)
call broadcast (remove_zero_projection)
call broadcast (cfac)
call broadcast (krook_odd)
if (.not.include_collisions) collisions_implicit = .false.
end subroutine read_parameters
subroutine init_collisions
use species, only: spec, nspec
use vpamu_grids, only: dvpa, dmu, mu, nmu
use stella_geometry, only: bmag
use stella_layouts
implicit none
integer :: is
real :: vnew_max
if (collisions_initialized) return
collisions_initialized = .true.
if (collision_model == "dougherty") then
write(*,*)
write(*,*) 'Coll. model: Dougherty'
if (collisions_implicit) then
write(*,*) 'Coll. algorithm: implicit'
if (vpa_operator) then
call init_vpadiff_matrix
call init_vpadiff_conserve
end if
if (mu_operator) then
call init_mudiff_matrix
call init_mudiff_conserve
end if
else
write(*,*) 'Coll. algorithm: explicit'
vnew_max = 0.0
do is = 1, nspec
vnew_max = max(vnew_max,maxval(spec(is)%vnew))
end do
cfl_dt_vpadiff = 2.0*dvpa**2/vnew_max
cfl_dt_mudiff = minval(bmag)/(vnew_max*maxval(mu(2:)/dmu(:nmu-1)**2))
end if
end if
if (collision_model == "fokker-planck") then
write(*,*) 'Coll. model: linearized Fokker-Planck'
write(*,*) 'Note: tested for linear CBC w. adiabatic e-'
call init_nusDpa
if (collisions_implicit) then
write(*,*) 'Coll. algorithm: implicit'
call init_fp_diffmatrix
call init_fp_conserve
else
write(*,*) 'Coll. algorithm: explicit'
vnew_max = 0.0
do is = 1, nspec
vnew_max = max(vnew_max,maxval(spec(is)%vnew))
end do
cfl_dt_vpadiff = 2.0*dvpa**2/vnew_max
cfl_dt_mudiff = minval(bmag)/(vnew_max*maxval(mu(2:)/dmu(:nmu-1)**2))
end if
end if
write(*,*)
end subroutine init_collisions
subroutine init_nusDpa
! AVB: compute the collision frequencies nuD, nus and nupa
use zgrid, only: nzgrid
use vpamu_grids, only: nmu, mu, dmu, vpa, nvpa
use stella_geometry, only: bmag
use species, only: spec
use spfunc, only: erf => erf_ext
use finite_differences, only: fd3pt
use vpamu_grids, only: maxwell_mu, maxwell_vpa
use constants, only: pi
implicit none
integer :: ia, imu, iv, iz, is
real :: x, Gf
if (.not.allocated(nus)) allocate (nus(nvpa,nmu,-nzgrid:nzgrid))
if (.not.allocated(nuD)) allocate (nuD(nvpa,nmu,-nzgrid:nzgrid))
if (.not.allocated(nupa)) allocate (nupa(nvpa,nmu,-nzgrid:nzgrid))
if (.not.allocated(nux)) allocate (nux(nvpa,nmu,-nzgrid:nzgrid))
if (.not.allocated(mw)) allocate (mw(nvpa,nmu,-nzgrid:nzgrid))
ia = 1
is = 1
do iz = -nzgrid, nzgrid
do iv = 1, nvpa
do imu = 1, nmu
x = sqrt(vpa(iv)**2 + 2*bmag(ia,iz)*mu(imu))
Gf = (erf(x)-x*(2/sqrt(pi))*exp(-x**2)) / (2*x**2)
nuD(iv,imu,iz) = spec(is)%vnew(is)*(erf(x)-Gf)/x**3
nus(iv,imu,iz) = spec(is)%vnew(is)*4*Gf/x
nupa(iv,imu,iz)= spec(is)%vnew(is)*2*Gf/x**3
mw(iv,imu,iz) = maxwell_vpa(iv,is)*maxwell_mu(1,iz,imu,is)
end do
end do
end do
nux = nupa - nuD
end subroutine init_nusDpa
subroutine finish_nusDpa
implicit none
if (allocated(nus)) deallocate (nus)
if (allocated(nuD)) deallocate (nuD)
if (allocated(nupa)) deallocate (nupa)
if (allocated(nux)) deallocate (nux)
if (allocated(mw)) deallocate (mw)
end subroutine finish_nusDpa
subroutine init_fp_diffmatrix
use stella_time, only: code_dt
use species, only: nspec, spec
use zgrid, only: nzgrid
use vpamu_grids, only: dvpa, vpa, nvpa, mu, nmu, maxwell_mu, maxwell_vpa, dmu, equally_spaced_mu_grid
use stella_layouts, only: kxkyz_lo
use stella_layouts, only: iky_idx, ikx_idx, iz_idx, is_idx
use stella_geometry, only: bmag
use dist_fn_arrays, only: kperp2
use constants, only: pi
implicit none
integer :: ikxkyz, iky, ikx, iz, is
integer :: imu, ia, iv
integer :: nc, nb, lldab
real :: vpap, vpam, vfac, mum, mup
real :: xpv, xmv, nupapv, nupamv, nuDpv, nuDmv, mwpv, mwmv, gam_mu, gam_mum, gam_mup
real :: mwm, mwp, nuDm, nuDp, nupam, nupap, xm, xp
if (.not.allocated(aa_blcs)) allocate (aa_blcs(nvpa,nmu,nmu,-nzgrid:nzgrid,nspec))
if (.not.allocated(bb_blcs)) allocate (bb_blcs(nvpa,nmu,nmu, kxkyz_lo%llim_proc:kxkyz_lo%ulim_alloc))
if (.not.allocated(cc_blcs)) allocate (cc_blcs(nvpa,nmu,nmu,-nzgrid:nzgrid,nspec))
if (.not.allocated(dd_vecs)) allocate (dd_vecs(nvpa,nmu,-nzgrid:nzgrid,nspec))
if (.not.allocated(cdiffmat_band)) allocate (cdiffmat_band(3*(nmu+1)+1,nmu*nvpa,kxkyz_lo%llim_proc:kxkyz_lo%ulim_alloc))
if (.not.allocated(blockmatrix)) allocate (blockmatrix(nvpa*nmu,nvpa*nmu,kxkyz_lo%llim_proc:kxkyz_lo%ulim_alloc))
if (.not.allocated(ipiv)) allocate (ipiv(nvpa*nmu,kxkyz_lo%llim_proc:kxkyz_lo%ulim_alloc))
ia = 1
vfac = 1 ! zero vpa-operator, in beta
aa_blcs = 0.
bb_blcs = 0.
cc_blcs = 0.
dd_vecs = 0.
! AVB: construct difference matrix stored in block form
! AVB: this matrix is nvpa*nmu x nvpa*nmu
! AVB: aa_blcs stores subdiagonal blocks, bb_blcs diagonal blocks and cc_blcs the superdiagonal blocks
! AVB: aa_blcs(1,:,:) and cc_blcs(nmu,:,:) are never used
do imu = 1, nmu
do iv = 1, nvpa
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)
is = is_idx(kxkyz_lo,ikxkyz)
if (iv == 1) then
vpap = 0.5*(vpa(iv)+vpa(iv+1))
mwpv = exp(-vpap**2)*maxwell_mu(1,iz,imu,is)
xpv = sqrt(vpap**2 + 2*bmag(ia,iz)*mu(imu))
nupapv = vfac*spec(is)%vnew(is)*2*(erf(xpv)-xpv*(2/sqrt(pi))*exp(-xpv**2)) / (2*xpv**2)/xpv**3
nuDpv = vfac*spec(is)%vnew(is)*(erf(xpv)-(erf(xpv)-xpv*(2/sqrt(pi))*exp(-xpv**2)) / (2*xpv**2))/xpv**3
if (imu == 1) then
! one-sided difference for mu-derivative at imu=1:
cc_blcs(iv,imu,imu+1,iz,is) = -vfac*code_dt*vpa(iv+1)*mu(imu)*nux(iv+1,imu,iz)*mw(iv+1,imu,iz)/mw(iv+1,imu+1,iz) / (2*dvpa) / dmu(imu)
cc_blcs(iv,imu,imu,iz,is) = -code_dt*0.5*(nupapv*vpap**2 + 2*nuDpv*bmag(ia,iz)*mu(imu))*mwpv / dvpa**2 / mw(iv+1,imu,iz) &
+vfac*code_dt*vpa(iv+1)*mu(imu)*nux(iv+1,imu,iz)*mw(iv+1,imu,iz)/mw(iv+1,imu,iz) / (2*dvpa) / dmu(imu)
bb_blcs(iv,imu,imu,ikxkyz) = 1 + code_dt*(0.5*(nupapv*vpap**2 + 2*nuDpv*bmag(ia,iz)*mu(imu))*mwpv) / dvpa**2 / mw(iv,imu,iz)
! mu operator
if (mu_operator) then
! use ghost cell at mu_{0} = 0, where term vanishes, so dmu(0) = mu(1).
mup = 0.5*(mu(imu)+mu(imu+1))
mwp = maxwell_vpa(iv,is)*exp(-2*bmag(ia,iz)*mup)
xp = sqrt(vpa(iv)**2 + 2*bmag(ia,iz)*mup)
nuDp = spec(is)%vnew(is)*(erf(xp)-(erf(xp)-xp*(2/sqrt(pi))*exp(-xp**2)) / (2*xp**2)) / xp**3
nupap = spec(is)%vnew(is)*2*(erf(xp)-xp*(2/sqrt(pi))*exp(-xp**2)) / (2*xp**2) / xp**3
gam_mu = 2*(nupa(iv,imu,iz)*mu(imu)**2+nuD(iv,imu,iz)*vpa(iv)**2/(2*bmag(ia,iz))*mu(imu))*mw(iv,imu,iz)
gam_mup = 2*(nupap*mup**2+nuDp*vpa(iv)**2/(2*bmag(ia,iz))*mup)*mwp
bb_blcs(iv,imu,imu,ikxkyz) = bb_blcs(iv,imu,imu,ikxkyz) - code_dt*(-gam_mu*1/dmu(imu) * dmu(imu)/2./mu(imu) &
+(-gam_mup*1/dmu(imu) + gam_mu*1/dmu(imu)) * mu(imu)/(dmu(imu)/2.)) / mw(iv,imu,iz) / (dmu(imu)/2.+mu(imu))
bb_blcs(iv,imu,imu+1,ikxkyz)= bb_blcs(iv,imu,imu+1,ikxkyz) - code_dt*(gam_mu * 1/dmu(imu) * dmu(imu)/2./mu(imu) &
+(gam_mup*1/dmu(imu) - gam_mu*1/dmu(imu)) * mu(imu)/(dmu(imu)/2.)) / mw(iv,imu+1,iz) / (dmu(imu)/2.+mu(imu))
! mixed derivative:
cc_blcs(iv,imu,imu,iz,is) = cc_blcs(iv,imu,imu,iz,is) - code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv+1,imu,iz)*(dmu(imu)/mu(imu)-mu(imu)/dmu(imu))) / (mu(imu)+dmu(imu)) / dvpa
cc_blcs(iv,imu,imu+1,iz,is) = cc_blcs(iv,imu,imu+1,iz,is) - code_dt*(vpa(iv)*mu(imu+1)*nux(iv,imu+1,iz)*mw(iv,imu+1,iz)*1/mw(iv+1,imu+1,iz)* mu(imu)/dmu(imu)) / (mu(imu)+dmu(imu)) / dvpa
bb_blcs(iv,imu,imu,ikxkyz) = bb_blcs(iv,imu,imu,ikxkyz) + code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv,imu,iz)*(dmu(imu)/mu(imu)-mu(imu)/dmu(imu))) / (mu(imu)+dmu(imu)) / dvpa
bb_blcs(iv,imu,imu+1,ikxkyz)= bb_blcs(iv,imu,imu+1,ikxkyz) + code_dt*(vpa(iv)*mu(imu+1)*nux(iv,imu+1,iz)*mw(iv,imu+1,iz)*1/mw(iv,imu+1,iz)* mu(imu)/dmu(imu)) / (mu(imu)+dmu(imu)) / dvpa
end if
else if (imu == nmu) then
! AVB: one-sided difference for mu-derivative at imu=nmu:
cc_blcs(iv,imu,imu-1,iz,is) = +vfac*code_dt*vpa(iv+1)*mu(imu)*nux(iv+1,imu,iz)*mw(iv+1,imu,iz)/mw(iv+1,imu-1,iz) / (2*dvpa) / dmu(nmu-1)
cc_blcs(iv,imu,imu,iz,is) = -code_dt*0.5*(nupapv*vpap**2 + 2*nuDpv*bmag(ia,iz)*mu(imu))*mwpv / dvpa**2 / mw(iv+1,imu,iz) &
-vfac*code_dt*vpa(iv+1)*mu(imu)*nux(iv+1,imu,iz)*mw(iv+1,imu,iz)/mw(iv+1,imu,iz) / (2*dvpa) / dmu(nmu-1)
bb_blcs(iv,imu,imu,ikxkyz) = 1 + code_dt*(0.5*(nupapv*vpap**2 + 2*nuDpv*bmag(ia,iz)*mu(imu))*mwpv) / dvpa**2 / mw(iv,imu,iz)
! mu operator
if (mu_operator) then
mum = 0.5*(mu(imu)+mu(imu-1))
mwm = maxwell_vpa(iv,is)*exp(-2*bmag(ia,iz)*mum)
xm = sqrt(vpa(iv)**2 + 2*bmag(ia,iz)*mum)
nuDm = spec(is)%vnew(is)*(erf(xm)-(erf(xm)-xm*(2/sqrt(pi))*exp(-xm**2)) / (2*xm**2)) / xm**3
nupam = spec(is)%vnew(is)*2*(erf(xm)-xm*(2/sqrt(pi))*exp(-xm**2)) / (2*xm**2) / xm**3
gam_mu = 2*(nupa(iv,imu,iz)*mu(imu)**2+nuD(iv,imu,iz)*vpa(iv)**2/(2*bmag(ia,iz))*mu(imu))*mw(iv,imu,iz)
gam_mum = 2*(nupam*mum**2+nuDm*vpa(iv)**2/(2*bmag(ia,iz))*mum)*mwm
bb_blcs(iv,imu,imu,ikxkyz) = bb_blcs(iv,imu,imu,ikxkyz) - code_dt*((gam_mu/dmu(imu-1) - gam_mum/dmu(imu-1))*dmu(imu-1)/(dmu(imu-1)/2.) &
+ (-gam_mu/dmu(imu-1))*dmu(imu-1)/2./dmu(imu-1)) / mw(iv,imu,iz) / (dmu(imu-1)/2.+dmu(imu-1))
bb_blcs(iv,imu,imu-1,ikxkyz)= bb_blcs(iv,imu,imu-1,ikxkyz) - code_dt*((-gam_mu/dmu(imu-1) + gam_mum * 1/dmu(imu-1))*dmu(imu-1)/(dmu(imu-1)/2.) &
+gam_mu * 1/dmu(imu-1)*dmu(imu-1)/2./dmu(imu-1)) / mw(iv,imu-1,iz) / (dmu(imu-1)/2.+dmu(imu-1))
! mixed derivative:
cc_blcs(iv,imu,imu,iz,is) = cc_blcs(iv,imu,imu,iz,is) - code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv+1,imu,iz)*(dmu(imu-1)/dmu(imu-1)-dmu(imu-1)/dmu(imu-1))) / (dmu(imu-1)+dmu(imu-1)) / (dvpa)
cc_blcs(iv,imu,imu-1,iz,is) = cc_blcs(iv,imu,imu-1,iz,is) + code_dt*(vpa(iv)*mu(imu-1)*nux(iv,imu-1,iz)*mw(iv,imu-1,iz)*1/mw(iv+1,imu-1,iz)* dmu(imu-1)/dmu(imu-1)) / (dmu(imu-1)+dmu(imu-1)) / (dvpa)
bb_blcs(iv,imu,imu,ikxkyz) = bb_blcs(iv,imu,imu,ikxkyz) + code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv,imu,iz)*(dmu(imu-1)/dmu(imu-1)-dmu(imu-1)/dmu(imu-1))) / (dmu(imu-1)+dmu(imu-1)) / (dvpa)
bb_blcs(iv,imu,imu-1,ikxkyz)= bb_blcs(iv,imu,imu-1,ikxkyz) - code_dt*(vpa(iv)*mu(imu-1)*nux(iv,imu-1,iz)*mw(iv,imu-1,iz)*1/mw(iv,imu-1,iz)* dmu(imu-1)/dmu(imu-1)) / (dmu(imu-1)+dmu(imu-1)) / (dvpa)
end if
else
cc_blcs(iv,imu,imu-1,iz,is) = +vfac*code_dt*vpa(iv+1)*mu(imu)*nux(iv+1,imu,iz)*mw(iv+1,imu,iz)/mw(iv+1,imu-1,iz)*dmu(imu)/dmu(imu-1) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
cc_blcs(iv,imu,imu+1,iz,is) = -vfac*code_dt*vpa(iv+1)*mu(imu)*nux(iv+1,imu,iz)*mw(iv+1,imu,iz)/mw(iv+1,imu+1,iz)*dmu(imu-1)/dmu(imu) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
cc_blcs(iv,imu,imu,iz,is) = -code_dt*0.5*(nupapv*vpap**2 + 2*nuDpv*bmag(ia,iz)*mu(imu))*mwpv / dvpa**2 / mw(iv+1,imu,iz) &
-vfac*code_dt*vpa(iv+1)*mu(imu)*nux(iv+1,imu,iz)*(dmu(imu)/dmu(imu-1) - dmu(imu-1)/dmu(imu)) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
bb_blcs(iv,imu,imu,ikxkyz) = 1 + code_dt*(0.5*(nupapv*vpap**2 + 2*nuDpv*bmag(ia,iz)*mu(imu))*mwpv) / dvpa**2 / mw(iv,imu,iz)
! mu operator
if (mu_operator) then
! quantities at mu_{i+1/2} and mu_{i-1/2}:
mup = 0.5*(mu(imu)+mu(imu+1))
mum = 0.5*(mu(imu)+mu(imu-1))
mwp = maxwell_vpa(iv,is)*exp(-2*bmag(ia,iz)*mup)
mwm = maxwell_vpa(iv,is)*exp(-2*bmag(ia,iz)*mum)
xp = sqrt(vpa(iv)**2 + 2*bmag(ia,iz)*mup)
xm = sqrt(vpa(iv)**2 + 2*bmag(ia,iz)*mum)
nuDp = spec(is)%vnew(is)*(erf(xp)-(erf(xp)-xp*(2/sqrt(pi))*exp(-xp**2)) / (2*xp**2)) / xp**3
nuDm = spec(is)%vnew(is)*(erf(xm)-(erf(xm)-xm*(2/sqrt(pi))*exp(-xm**2)) / (2*xm**2)) / xm**3
nupap = spec(is)%vnew(is)*2*(erf(xp)-xp*(2/sqrt(pi))*exp(-xp**2)) / (2*xp**2) / xp**3
nupam = spec(is)%vnew(is)*2*(erf(xm)-xm*(2/sqrt(pi))*exp(-xm**2)) / (2*xm**2) / xm**3
gam_mu = 2*(nupa(iv,imu,iz)*mu(imu)**2+nuD(iv,imu,iz)*vpa(iv)**2/(2*bmag(ia,iz))*mu(imu))*mw(iv,imu,iz)
gam_mum = 2*(nupam*mum**2+nuDm*vpa(iv)**2/(2*bmag(ia,iz))*mum)*mwm
gam_mup = 2*(nupap*mup**2+nuDp*vpa(iv)**2/(2*bmag(ia,iz))*mup)*mwp
bb_blcs(iv,imu,imu,ikxkyz) = bb_blcs(iv,imu,imu,ikxkyz) - code_dt*((gam_mu*(dmu(imu)/dmu(imu-1) - dmu(imu-1)/dmu(imu)) / (dmu(imu-1)+dmu(imu)) - gam_mum / dmu(imu-1))*dmu(imu)/dmu(imu-1) &
+(-gam_mu*(dmu(imu)/dmu(imu-1) - dmu(imu-1)/dmu(imu)) / (dmu(imu-1)+dmu(imu)) - gam_mup / dmu(imu))* dmu(imu-1)/dmu(imu)) / mw(iv,imu,iz) * 2./(dmu(imu-1)+dmu(imu))
bb_blcs(iv,imu,imu-1,ikxkyz) = bb_blcs(iv,imu,imu-1,ikxkyz) - code_dt*((-gam_mu * 1*dmu(imu)/dmu(imu-1) / (dmu(imu-1)+dmu(imu)) + gam_mum * 1/dmu(imu-1)) * dmu(imu)/dmu(imu-1) &
+gam_mu * 1*dmu(imu)/dmu(imu-1) / (dmu(imu-1)+dmu(imu)) * dmu(imu-1)/dmu(imu)) / mw(iv,imu-1,iz) * 2./(dmu(imu-1)+dmu(imu))
bb_blcs(iv,imu,imu+1,ikxkyz) = bb_blcs(iv,imu,imu+1,ikxkyz) - code_dt*( gam_mu*dmu(imu-1)/dmu(imu) / (dmu(imu-1)+dmu(imu)) * dmu(imu)/dmu(imu-1) &
+ (gam_mup/dmu(imu) - gam_mu*dmu(imu-1)/dmu(imu) / (dmu(imu-1)+dmu(imu))) * dmu(imu-1)/dmu(imu)) / mw(iv,imu+1,iz) * 2/(dmu(imu-1)+dmu(imu))
! mixed derivative:
cc_blcs(iv,imu,imu,iz,is) = cc_blcs(iv,imu,imu,iz,is) - code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv+1,imu,iz)*(dmu(imu)/dmu(imu-1)-dmu(imu-1)/dmu(imu))) / (dmu(imu-1)+dmu(imu)) / (dvpa)
cc_blcs(iv,imu,imu-1,iz,is) = cc_blcs(iv,imu,imu-1,iz,is) + code_dt*(vpa(iv)*mu(imu-1)*nux(iv,imu-1,iz)*mw(iv,imu-1,iz)*1/mw(iv+1,imu-1,iz)* dmu(imu)/dmu(imu-1)) / (dmu(imu-1)+dmu(imu)) / (dvpa)
cc_blcs(iv,imu,imu+1,iz,is) = cc_blcs(iv,imu,imu+1,iz,is) - code_dt*(vpa(iv)*mu(imu+1)*nux(iv,imu+1,iz)*mw(iv,imu+1,iz)*1/mw(iv+1,imu+1,iz)* dmu(imu-1)/dmu(imu)) / (dmu(imu-1)+dmu(imu)) / (dvpa)
bb_blcs(iv,imu,imu,ikxkyz) = bb_blcs(iv,imu,imu,ikxkyz) + code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv,imu,iz)*(dmu(imu)/dmu(imu-1)-dmu(imu-1)/dmu(imu))) / (dmu(imu-1)+dmu(imu)) / (dvpa)
bb_blcs(iv,imu,imu-1,ikxkyz) = bb_blcs(iv,imu,imu-1,ikxkyz) - code_dt*(vpa(iv)*mu(imu-1)*nux(iv,imu-1,iz)*mw(iv,imu-1,iz)*1/mw(iv,imu-1,iz)* dmu(imu)/dmu(imu-1)) / (dmu(imu-1)+dmu(imu)) / (dvpa)
bb_blcs(iv,imu,imu+1,ikxkyz) = bb_blcs(iv,imu,imu+1,ikxkyz) + code_dt*(vpa(iv)*mu(imu+1)*nux(iv,imu+1,iz)*mw(iv,imu+1,iz)*1/mw(iv,imu+1,iz)* dmu(imu-1)/dmu(imu)) / (dmu(imu-1)+dmu(imu)) / (dvpa)
end if
end if
else if (iv == nvpa) then
vpam = 0.5*(vpa(iv)+vpa(iv-1))
mwmv = exp(-vpam**2)*maxwell_mu(1,iz,imu,is)
xmv = sqrt(vpam**2 + 2*bmag(ia,iz)*mu(imu))
nupamv = vfac*spec(is)%vnew(is)*2*(erf(xmv)-xmv*(2/sqrt(pi))*exp(-xmv**2)) / (2*xmv**2)/xmv**3
nuDmv = vfac*spec(is)%vnew(is)*(erf(xmv)-(erf(xmv)-xmv*(2/sqrt(pi))*exp(-xmv**2)) / (2*xmv**2))/xmv**3
if (imu == 1) then
! one-sided difference for mu-derivative at imu=1:
aa_blcs(iv,imu,imu+1,iz,is) = +vfac*code_dt*vpa(iv-1)*mu(imu)*nux(iv-1,imu,iz)*mw(iv-1,imu,iz)/mw(iv-1,imu+1,iz) / (2*dvpa) / dmu(imu)
aa_blcs(iv,imu,imu,iz,is) = -code_dt*0.5*(nupamv*vpam**2 + 2*nuDmv*bmag(ia,iz)*mu(imu))*mwmv / dvpa**2 / mw(iv-1,imu,iz) &
-vfac*code_dt*vpa(iv-1)*mu(imu)*nux(iv-1,imu,iz)*mw(iv-1,imu,iz)/mw(iv-1,imu,iz) / (2*dvpa) / dmu(imu)
bb_blcs(iv,imu,imu,ikxkyz) = 1 + code_dt*(0.5*(nupamv*vpam**2 + 2*nuDmv*bmag(ia,iz)*mu(imu))*mwmv) / dvpa**2 / mw(iv,imu,iz)
! mu operator
if (mu_operator) then
! use ghost cell at mu_{0} = 0, where term vanishes, so dmu(0) = mu(1).
mup = 0.5*(mu(imu)+mu(imu+1))
mwp = maxwell_vpa(iv,is)*exp(-2*bmag(ia,iz)*mup)
xp = sqrt(vpa(iv)**2 + 2*bmag(ia,iz)*mup)
nuDp = spec(is)%vnew(is)*(erf(xp)-(erf(xp)-xp*(2/sqrt(pi))*exp(-xp**2)) / (2*xp**2)) / xp**3
nupap = spec(is)%vnew(is)*2*(erf(xp)-xp*(2/sqrt(pi))*exp(-xp**2)) / (2*xp**2) / xp**3
gam_mu = 2*(nupa(iv,imu,iz)*mu(imu)**2+nuD(iv,imu,iz)*vpa(iv)**2/(2*bmag(ia,iz))*mu(imu))*mw(iv,imu,iz)
gam_mup = 2*(nupap*mup**2+nuDp*vpa(iv)**2/(2*bmag(ia,iz))*mup)*mwp
bb_blcs(iv,imu,imu,ikxkyz) = bb_blcs(iv,imu,imu,ikxkyz) - code_dt*(-gam_mu *1/dmu(imu) * dmu(imu)/2./mu(imu) &
+(-gam_mup*1/dmu(imu) + gam_mu*1/dmu(imu)) * mu(imu)/(dmu(imu)/2.)) / mw(iv,imu,iz) / (dmu(imu)/2.+mu(imu))
bb_blcs(iv,imu,imu+1,ikxkyz) = bb_blcs(iv,imu,imu+1,ikxkyz) - code_dt*(gam_mu * 1/dmu(imu) * dmu(imu)/2./mu(imu) &
+(gam_mup* 1/dmu(imu) - gam_mu*1/dmu(imu)) * mu(imu)/(dmu(imu)/2.)) / mw(iv,imu+1,iz) / (dmu(imu)/2.+mu(imu))
! mixed derivative:
aa_blcs(iv,imu,imu,iz,is) = aa_blcs(iv,imu,imu,iz,is) + code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv-1,imu,iz)*(dmu(imu)/mu(imu)-mu(imu)/dmu(imu))) / (mu(imu)+dmu(imu)) / (dvpa)
aa_blcs(iv,imu,imu+1,iz,is) = aa_blcs(iv,imu,imu+1,iz,is) + code_dt*(vpa(iv)*mu(imu+1)*nux(iv,imu+1,iz)*mw(iv,imu+1,iz)*1/mw(iv-1,imu+1,iz)* mu(imu)/dmu(imu)) / (mu(imu)+dmu(imu)) / (dvpa)
bb_blcs(iv,imu,imu,ikxkyz) = bb_blcs(iv,imu,imu,ikxkyz) - code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv,imu,iz)*(dmu(imu)/mu(imu)-mu(imu)/dmu(imu))) / (mu(imu)+dmu(imu)) / (dvpa)
bb_blcs(iv,imu,imu+1,ikxkyz) = bb_blcs(iv,imu,imu+1,ikxkyz) - code_dt*(vpa(iv)*mu(imu+1)*nux(iv,imu+1,iz)*mw(iv,imu+1,iz)*1/mw(iv,imu+1,iz)* mu(imu)/dmu(imu)) / (mu(imu)+dmu(imu)) / (dvpa)
end if
else if (imu == nmu) then
! one-sided difference for mu-derivative at imu=nmu:
aa_blcs(iv,imu,imu-1,iz,is) = -vfac*code_dt*vpa(iv-1)*mu(imu)*nux(iv-1,imu,iz)*mw(iv-1,imu,iz)/mw(iv-1,imu-1,iz) / (2*dvpa) / dmu(nmu-1)
aa_blcs(iv,imu,imu,iz,is) = -code_dt*0.5*(nupamv*vpam**2 + 2*nuDmv*bmag(ia,iz)*mu(imu))*mwmv / dvpa**2 / mw(iv-1,imu,iz) &
+vfac*code_dt*vpa(iv-1)*mu(imu)*nux(iv-1,imu,iz)*mw(iv-1,imu,iz)/mw(iv-1,imu,iz) / (2*dvpa) / dmu(nmu-1)
bb_blcs(iv,imu,imu,ikxkyz) = 1 + code_dt*(0.5*(nupamv*vpam**2 + 2*nuDmv*bmag(ia,iz)*mu(imu))*mwmv) / dvpa**2 / mw(iv,imu,iz)
! mu operator
if (mu_operator) then
mum = 0.5*(mu(imu)+mu(imu-1))
mwm = maxwell_vpa(iv,is)*exp(-2*bmag(ia,iz)*mum)
xm = sqrt(vpa(iv)**2 + 2*bmag(ia,iz)*mum)
nuDm = spec(is)%vnew(is)*(erf(xm)-(erf(xm)-xm*(2/sqrt(pi))*exp(-xm**2)) / (2*xm**2)) / xm**3
nupam = spec(is)%vnew(is)*2*(erf(xm)-xm*(2/sqrt(pi))*exp(-xm**2)) / (2*xm**2) / xm**3
gam_mu = 2*(nupa(iv,imu,iz)*mu(imu)**2+nuD(iv,imu,iz)*vpa(iv)**2/(2*bmag(ia,iz))*mu(imu))*mw(iv,imu,iz)
gam_mum = 2*(nupam*mum**2+nuDm*vpa(iv)**2/(2*bmag(ia,iz))*mum)*mwm
bb_blcs(iv,imu,imu,ikxkyz) = bb_blcs(iv,imu,imu,ikxkyz) - code_dt*((gam_mu/dmu(imu-1) - gam_mum/dmu(imu-1))*dmu(imu-1)/(dmu(imu-1)/2.) &
+ (-gam_mu/dmu(imu-1))*dmu(imu-1)/2./dmu(imu-1)) / mw(iv,imu,iz) / (dmu(imu-1)/2.+dmu(imu-1))
bb_blcs(iv,imu,imu-1,ikxkyz) = bb_blcs(iv,imu,imu-1,ikxkyz) - code_dt*((-gam_mu/dmu(imu-1) + gam_mum * 1/dmu(imu-1))*dmu(imu-1)/(dmu(imu-1)/2.) &
+gam_mu * 1/dmu(imu-1)*dmu(imu-1)/2./dmu(imu-1)) / mw(iv,imu-1,iz) / (dmu(imu-1)/2.+dmu(imu-1))
! mixed derivative:
aa_blcs(iv,imu,imu,iz,is) = aa_blcs(iv,imu,imu,iz,is) + code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv-1,imu,iz)*(dmu(imu-1)/dmu(imu-1)-dmu(imu-1)/dmu(imu-1))) / (dmu(imu-1)+dmu(imu-1)) / (dvpa)
aa_blcs(iv,imu,imu-1,iz,is) = aa_blcs(iv,imu,imu-1,iz,is) - code_dt*(vpa(iv)*mu(imu-1)*nux(iv,imu-1,iz)*mw(iv,imu-1,iz)*1/mw(iv-1,imu-1,iz)* dmu(imu-1)/dmu(imu-1)) / (dmu(imu-1)+dmu(imu-1)) / (dvpa)
bb_blcs(iv,imu,imu,ikxkyz) = bb_blcs(iv,imu,imu,ikxkyz) - code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv,imu,iz)*(dmu(imu-1)/dmu(imu-1)-dmu(imu-1)/dmu(imu-1))) / (dmu(imu-1)+dmu(imu-1)) / (dvpa)
bb_blcs(iv,imu,imu-1,ikxkyz) = bb_blcs(iv,imu,imu-1,ikxkyz) + code_dt*(vpa(iv)*mu(imu-1)*nux(iv,imu-1,iz)*mw(iv,imu-1,iz)*1/mw(iv,imu-1,iz)* dmu(imu-1)/dmu(imu-1)) / (dmu(imu-1)+dmu(imu-1)) / (dvpa)
end if
else
aa_blcs(iv,imu,imu-1,iz,is) = -vfac*code_dt*vpa(iv-1)*mu(imu)*nux(iv-1,imu,iz)*mw(iv-1,imu,iz)/mw(iv-1,imu-1,iz)* dmu(imu)/dmu(imu-1) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
aa_blcs(iv,imu,imu+1,iz,is) = +vfac*code_dt*vpa(iv-1)*mu(imu)*nux(iv-1,imu,iz)*mw(iv-1,imu,iz)/mw(iv-1,imu+1,iz)* dmu(imu-1)/dmu(imu) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
aa_blcs(iv,imu,imu,iz,is) = -code_dt*0.5*(nupamv*vpam**2 + 2*nuDmv*bmag(ia,iz)*mu(imu))*mwmv / dvpa**2 / mw(iv-1,imu,iz) &
+vfac*code_dt*vpa(iv-1)*mu(imu)*nux(iv-1,imu,iz)*(dmu(imu)/dmu(imu-1) - dmu(imu-1)/dmu(imu)) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
bb_blcs(iv,imu,imu,ikxkyz) = 1 + code_dt*(0.5*(nupamv*vpam**2 + 2*nuDmv*bmag(ia,iz)*mu(imu))*mwmv) / dvpa**2 / mw(iv,imu,iz)
! mu operator
if (mu_operator) then
! quantities at mu_{i+1/2} and mu_{i-1/2}:
mup = 0.5*(mu(imu)+mu(imu+1))
mum = 0.5*(mu(imu)+mu(imu-1))
mwp = maxwell_vpa(iv,is)*exp(-2*bmag(ia,iz)*mup)
mwm = maxwell_vpa(iv,is)*exp(-2*bmag(ia,iz)*mum)
xp = sqrt(vpa(iv)**2 + 2*bmag(ia,iz)*mup)
xm = sqrt(vpa(iv)**2 + 2*bmag(ia,iz)*mum)
nuDp = spec(is)%vnew(is)*(erf(xp)-(erf(xp)-xp*(2/sqrt(pi))*exp(-xp**2)) / (2*xp**2)) / xp**3
nuDm = spec(is)%vnew(is)*(erf(xm)-(erf(xm)-xm*(2/sqrt(pi))*exp(-xm**2)) / (2*xm**2)) / xm**3
nupap = spec(is)%vnew(is)*2*(erf(xp)-xp*(2/sqrt(pi))*exp(-xp**2)) / (2*xp**2) / xp**3
nupam = spec(is)%vnew(is)*2*(erf(xm)-xm*(2/sqrt(pi))*exp(-xm**2)) / (2*xm**2) / xm**3
gam_mu = 2*(nupa(iv,imu,iz)*mu(imu)**2+nuD(iv,imu,iz)*vpa(iv)**2/(2*bmag(ia,iz))*mu(imu))*mw(iv,imu,iz)
gam_mum = 2*(nupam*mum**2+nuDm*vpa(iv)**2/(2*bmag(ia,iz))*mum)*mwm
gam_mup = 2*(nupap*mup**2+nuDp*vpa(iv)**2/(2*bmag(ia,iz))*mup)*mwp
bb_blcs(iv,imu,imu,ikxkyz) = bb_blcs(iv,imu,imu,ikxkyz) - code_dt*((gam_mu*(dmu(imu)/dmu(imu-1) - dmu(imu-1)/dmu(imu)) / (dmu(imu-1)+dmu(imu)) - gam_mum / dmu(imu-1))*dmu(imu)/dmu(imu-1) &
+(-gam_mu*(dmu(imu)/dmu(imu-1) - dmu(imu-1)/dmu(imu)) / (dmu(imu-1)+dmu(imu)) - gam_mup / dmu(imu))* dmu(imu-1)/dmu(imu)) / mw(iv,imu,iz) * 2./(dmu(imu-1)+dmu(imu))
bb_blcs(iv,imu,imu-1,ikxkyz) = bb_blcs(iv,imu,imu-1,ikxkyz) - code_dt*((-gam_mu * 1*dmu(imu)/dmu(imu-1) / (dmu(imu-1)+dmu(imu)) + gam_mum * 1/dmu(imu-1)) * dmu(imu)/dmu(imu-1) &
+gam_mu * 1*dmu(imu)/dmu(imu-1) / (dmu(imu-1)+dmu(imu)) * dmu(imu-1)/dmu(imu)) / mw(iv,imu-1,iz) * 2./(dmu(imu-1)+dmu(imu))
bb_blcs(iv,imu,imu+1,ikxkyz) = bb_blcs(iv,imu,imu+1,ikxkyz) - code_dt*( gam_mu*dmu(imu-1)/dmu(imu) / (dmu(imu-1)+dmu(imu)) * dmu(imu)/dmu(imu-1) &
+ (gam_mup/dmu(imu) - gam_mu*dmu(imu-1)/dmu(imu) / (dmu(imu-1)+dmu(imu))) * dmu(imu-1)/dmu(imu)) / mw(iv,imu+1,iz) * 2/(dmu(imu-1)+dmu(imu))
! mixed derivative, one-sided difference in vpa at iv = nvpa:
aa_blcs(iv,imu,imu,iz,is) = aa_blcs(iv,imu,imu,iz,is) + code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv-1,imu,iz)*(dmu(imu)/dmu(imu-1)-dmu(imu-1)/dmu(imu))) / (dmu(imu-1)+dmu(imu)) / (dvpa)
aa_blcs(iv,imu,imu-1,iz,is) = aa_blcs(iv,imu,imu-1,iz,is) - code_dt*(vpa(iv)*mu(imu-1)*nux(iv,imu-1,iz)*mw(iv,imu-1,iz)*1/mw(iv-1,imu-1,iz)* dmu(imu)/dmu(imu-1)) / (dmu(imu-1)+dmu(imu)) / (dvpa)
aa_blcs(iv,imu,imu+1,iz,is) = aa_blcs(iv,imu,imu+1,iz,is) + code_dt*(vpa(iv)*mu(imu+1)*nux(iv,imu+1,iz)*mw(iv,imu+1,iz)*1/mw(iv-1,imu+1,iz)* dmu(imu-1)/dmu(imu)) / (dmu(imu-1)+dmu(imu)) / (dvpa)
bb_blcs(iv,imu,imu,ikxkyz) = bb_blcs(iv,imu,imu,ikxkyz) - code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv,imu,iz)*(dmu(imu)/dmu(imu-1)-dmu(imu-1)/dmu(imu))) / (dmu(imu-1)+dmu(imu)) / (dvpa)
bb_blcs(iv,imu,imu-1,ikxkyz) = bb_blcs(iv,imu,imu-1,ikxkyz) + code_dt*(vpa(iv)*mu(imu-1)*nux(iv,imu-1,iz)*mw(iv,imu-1,iz)*1/mw(iv,imu-1,iz)* dmu(imu)/dmu(imu-1)) / (dmu(imu-1)+dmu(imu)) / (dvpa)
bb_blcs(iv,imu,imu+1,ikxkyz) = bb_blcs(iv,imu,imu+1,ikxkyz) - code_dt*(vpa(iv)*mu(imu+1)*nux(iv,imu+1,iz)*mw(iv,imu+1,iz)*1/mw(iv,imu+1,iz)* dmu(imu-1)/dmu(imu)) / (dmu(imu-1)+dmu(imu)) / (dvpa)
end if
end if
else
vpam = 0.5*(vpa(iv)+vpa(iv-1))
vpap = 0.5*(vpa(iv)+vpa(iv+1))
mwmv = exp(-vpam**2)*maxwell_mu(1,iz,imu,is)
mwpv = exp(-vpap**2)*maxwell_mu(1,iz,imu,is)
xpv = sqrt(vpap**2 + 2*bmag(ia,iz)*mu(imu))
xmv = sqrt(vpam**2 + 2*bmag(ia,iz)*mu(imu))
nupamv = vfac*spec(is)%vnew(is)*2*(erf(xmv)-xmv*(2/sqrt(pi))*exp(-xmv**2)) / (2*xmv**2)/xmv**3
nupapv = vfac*spec(is)%vnew(is)*2*(erf(xpv)-xpv*(2/sqrt(pi))*exp(-xpv**2)) / (2*xpv**2)/xpv**3
nuDmv = vfac*spec(is)%vnew(is)*(erf(xmv)-(erf(xmv)-xmv*(2/sqrt(pi))*exp(-xmv**2)) / (2*xmv**2))/xmv**3
nuDpv = vfac*spec(is)%vnew(is)*(erf(xpv)-(erf(xpv)-xpv*(2/sqrt(pi))*exp(-xpv**2)) / (2*xpv**2))/xpv**3
if (imu == 1) then
! one-sided difference for mu-derivative at imu=1:
aa_blcs(iv,imu,imu+1,iz,is) = +vfac*code_dt*vpa(iv-1)*mu(imu)*nux(iv-1,imu,iz)*mw(iv-1,imu,iz)/mw(iv-1,imu+1,iz) / (2*dvpa) / dmu(imu)
aa_blcs(iv,imu,imu,iz,is) = -code_dt*0.5*(nupamv*vpam**2 + 2*nuDmv*bmag(ia,iz)*mu(imu))*mwmv / dvpa**2 / mw(iv-1,imu,iz) -vfac*code_dt*vpa(iv-1)*mu(imu)*nux(iv-1,imu,iz) / (2*dvpa) / dmu(imu)
cc_blcs(iv,imu,imu+1,iz,is) = -vfac*code_dt*vpa(iv+1)*mu(imu)*nux(iv+1,imu,iz)*mw(iv+1,imu,iz)/mw(iv+1,imu+1,iz) / (2*dvpa) / dmu(imu)
cc_blcs(iv,imu,imu,iz,is) = -code_dt*0.5*(nupapv*vpap**2 + 2*nuDpv*bmag(ia,iz)*mu(imu))*mwpv / dvpa**2 / mw(iv+1,imu,iz) +vfac*code_dt*vpa(iv+1)*mu(imu)*nux(iv+1,imu,iz) / (2*dvpa) / dmu(imu)
bb_blcs(iv,imu,imu,ikxkyz) = 1 + code_dt*(0.5*(nupapv*vpap**2 + 2*nuDpv*bmag(ia,iz)*mu(imu))*mwpv + 0.5*(nupamv*vpam**2 + 2*nuDmv*bmag(ia,iz)*mu(imu))*mwmv) / dvpa**2 / mw(iv,imu,iz)
! mu operator
if (mu_operator) then
! use ghost cell at mu_{0} = 0, where term vanishes, so dmu(0) = mu(1).
mup = 0.5*(mu(imu)+mu(imu+1))
mwp = maxwell_vpa(iv,is)*exp(-2*bmag(ia,iz)*mup)
xp = sqrt(vpa(iv)**2 + 2*bmag(ia,iz)*mup)
nuDp = spec(is)%vnew(is)*(erf(xp)-(erf(xp)-xp*(2/sqrt(pi))*exp(-xp**2)) / (2*xp**2)) / xp**3
nupap = spec(is)%vnew(is)*2*(erf(xp)-xp*(2/sqrt(pi))*exp(-xp**2)) / (2*xp**2) / xp**3
gam_mu = 2*(nupa(iv,imu,iz)*mu(imu)**2+nuD(iv,imu,iz)*vpa(iv)**2/(2*bmag(ia,iz))*mu(imu))*mw(iv,imu,iz)
gam_mup = 2*(nupap*mup**2+nuDp*vpa(iv)**2/(2*bmag(ia,iz))*mup)*mwp
bb_blcs(iv,imu,imu,ikxkyz) = bb_blcs(iv,imu,imu,ikxkyz) - code_dt*(-gam_mu *1/dmu(imu) * dmu(imu)/2./mu(imu) &
+(-gam_mup*1/dmu(imu) + gam_mu*1/dmu(imu)) * mu(imu)/(dmu(imu)/2.)) / mw(iv,imu,iz) / (dmu(imu)/2.+mu(imu))
bb_blcs(iv,imu,imu+1,ikxkyz) = bb_blcs(iv,imu,imu+1,ikxkyz) - code_dt*(gam_mu * 1/dmu(imu) * dmu(imu)/2./mu(imu) &
+(gam_mup* 1/dmu(imu) - gam_mu*1/dmu(imu)) * mu(imu)/(dmu(imu)/2.)) / mw(iv,imu+1,iz) / (dmu(imu)/2.+mu(imu))
! mixed derivative:
aa_blcs(iv,imu,imu,iz,is) = aa_blcs(iv,imu,imu,iz,is) + code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv-1,imu,iz)*(dmu(imu)/mu(imu)-mu(imu)/dmu(imu))) / (mu(imu)+dmu(imu)) / (2*dvpa)
aa_blcs(iv,imu,imu+1,iz,is) = aa_blcs(iv,imu,imu+1,iz,is) + code_dt*(vpa(iv)*mu(imu+1)*nux(iv,imu+1,iz)*mw(iv,imu+1,iz)*1/mw(iv-1,imu+1,iz)* mu(imu)/dmu(imu)) / (mu(imu)+dmu(imu)) / (2*dvpa)
cc_blcs(iv,imu,imu,iz,is) = cc_blcs(iv,imu,imu,iz,is) - code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv+1,imu,iz)*(dmu(imu)/mu(imu)-mu(imu)/dmu(imu))) / (mu(imu)+dmu(imu)) / (2*dvpa)
cc_blcs(iv,imu,imu+1,iz,is) = cc_blcs(iv,imu,imu+1,iz,is) - code_dt*(vpa(iv)*mu(imu+1)*nux(iv,imu+1,iz)*mw(iv,imu+1,iz)*1/mw(iv+1,imu+1,iz)* mu(imu)/dmu(imu)) / (mu(imu)+dmu(imu)) / (2*dvpa)
end if
else if (imu == nmu) then
! one-sided difference for mu-derivative at imu=nmu:
aa_blcs(iv,imu,imu-1,iz,is) = -vfac*code_dt*vpa(iv-1)*mu(imu)*nux(iv-1,imu,iz)*mw(iv-1,imu,iz)/mw(iv-1,imu-1,iz) / (2*dvpa) / dmu(nmu-1)
aa_blcs(iv,imu,imu,iz,is) = -code_dt*0.5*(nupamv*vpam**2 + 2*nuDmv*bmag(ia,iz)*mu(imu))*mwmv / dvpa**2 / mw(iv-1,imu,iz) + vfac*code_dt*vpa(iv-1)*mu(imu)*nux(iv-1,imu,iz) / (2*dvpa) / dmu(nmu-1)
cc_blcs(iv,imu,imu-1,iz,is) = +vfac*code_dt*vpa(iv+1)*mu(imu)*nux(iv+1,imu,iz)*mw(iv+1,imu,iz)/mw(iv+1,imu-1,iz) / (2*dvpa) / dmu(nmu-1)
cc_blcs(iv,imu,imu,iz,is) = -code_dt*0.5*(nupapv*vpap**2 + 2*nuDpv*bmag(ia,iz)*mu(imu))*mwpv / dvpa**2 / mw(iv+1,imu,iz) - vfac*code_dt*vpa(iv+1)*mu(imu)*nux(iv+1,imu,iz) / (2*dvpa) / dmu(nmu-1)
bb_blcs(iv,imu,imu,ikxkyz) = 1 + code_dt*(0.5*(nupapv*vpap**2 + 2*nuDpv*bmag(ia,iz)*mu(imu))*mwpv + 0.5*(nupamv*vpam**2 + 2*nuDmv*bmag(ia,iz)*mu(imu))*mwmv) / dvpa**2 / mw(iv,imu,iz)
! mu operator
if (mu_operator) then
mum = 0.5*(mu(imu)+mu(imu-1))
mwm = maxwell_vpa(iv,is)*exp(-2*bmag(ia,iz)*mum)
xm = sqrt(vpa(iv)**2 + 2*bmag(ia,iz)*mum)
nuDm = spec(is)%vnew(is)*(erf(xm)-(erf(xm)-xm*(2/sqrt(pi))*exp(-xm**2)) / (2*xm**2)) / xm**3
nupam = spec(is)%vnew(is)*2*(erf(xm)-xm*(2/sqrt(pi))*exp(-xm**2)) / (2*xm**2) / xm**3
gam_mu = 2*(nupa(iv,imu,iz)*mu(imu)**2+nuD(iv,imu,iz)*vpa(iv)**2/(2*bmag(ia,iz))*mu(imu))*mw(iv,imu,iz)
gam_mum = 2*(nupam*mum**2+nuDm*vpa(iv)**2/(2*bmag(ia,iz))*mum)*mwm
bb_blcs(iv,imu,imu,ikxkyz) = bb_blcs(iv,imu,imu,ikxkyz) - code_dt*((gam_mu/dmu(imu-1) - gam_mum/dmu(imu-1))*dmu(imu-1)/(dmu(imu-1)/2.) &
+ (-gam_mu/dmu(imu-1))*dmu(imu-1)/2./dmu(imu-1)) / mw(iv,imu,iz) / (dmu(imu-1)/2.+dmu(imu-1))
bb_blcs(iv,imu,imu-1,ikxkyz) = bb_blcs(iv,imu,imu-1,ikxkyz) - code_dt*((-gam_mu/dmu(imu-1) + gam_mum * 1/dmu(imu-1))*dmu(imu-1)/(dmu(imu-1)/2.) &
+gam_mu * 1/dmu(imu-1)*dmu(imu-1)/2./dmu(imu-1)) / mw(iv,imu-1,iz) / (dmu(imu-1)/2.+dmu(imu-1))
! mixed derivative:
aa_blcs(iv,imu,imu,iz,is) = aa_blcs(iv,imu,imu,iz,is) + code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv-1,imu,iz)*(dmu(imu-1)/dmu(imu-1)-dmu(imu-1)/dmu(imu-1))) / (dmu(imu-1)+dmu(imu-1)) / (2*dvpa)
aa_blcs(iv,imu,imu-1,iz,is) = aa_blcs(iv,imu,imu-1,iz,is) - code_dt*(vpa(iv)*mu(imu-1)*nux(iv,imu-1,iz)*mw(iv,imu-1,iz)*1/mw(iv-1,imu-1,iz)* dmu(imu-1)/dmu(imu-1)) / (dmu(imu-1)+dmu(imu-1)) / (2*dvpa)
cc_blcs(iv,imu,imu,iz,is) = cc_blcs(iv,imu,imu,iz,is) - code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv+1,imu,iz)*(dmu(imu-1)/dmu(imu-1)-dmu(imu-1)/dmu(imu-1))) / (dmu(imu-1)+dmu(imu-1)) / (2*dvpa)
cc_blcs(iv,imu,imu-1,iz,is) = cc_blcs(iv,imu,imu-1,iz,is) + code_dt*(vpa(iv)*mu(imu-1)*nux(iv,imu-1,iz)*mw(iv,imu-1,iz)*1/mw(iv+1,imu-1,iz)* dmu(imu-1)/dmu(imu-1)) / (dmu(imu-1)+dmu(imu-1)) / (2*dvpa)
end if
else
! vpa operator (interior treatment):
aa_blcs(iv,imu,imu,iz,is) = -code_dt*0.5*(nupamv*vpam**2 + 2*nuDmv*bmag(ia,iz)*mu(imu))*mwmv / dvpa**2 / mw(iv-1,imu,iz) &
+vfac*code_dt*vpa(iv-1)*mu(imu)*nux(iv-1,imu,iz) * (dmu(imu)/dmu(imu-1) - dmu(imu-1)/dmu(imu)) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
cc_blcs(iv,imu,imu,iz,is) = -code_dt*0.5*(nupapv*vpap**2 + 2*nuDpv*bmag(ia,iz)*mu(imu))*mwpv / dvpa**2 / mw(iv+1,imu,iz) &
-vfac*code_dt*vpa(iv+1)*mu(imu)*nux(iv+1,imu,iz) * (dmu(imu)/dmu(imu-1) - dmu(imu-1)/dmu(imu)) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
bb_blcs(iv,imu,imu,ikxkyz) = 1 + code_dt*(0.5*(nupapv*vpap**2 + 2*nuDpv*bmag(ia,iz)*mu(imu))*mwpv + 0.5*(nupamv*vpam**2 + 2*nuDmv*bmag(ia,iz)*mu(imu))*mwmv) / dvpa**2 / mw(iv,imu,iz)
! vpa operator, mixed (interior treatment):
aa_blcs(iv,imu,imu-1,iz,is) = -vfac*code_dt*vpa(iv-1)*mu(imu)*nux(iv-1,imu,iz)*mw(iv-1,imu,iz)/mw(iv-1,imu-1,iz)* dmu(imu)/dmu(imu-1) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
aa_blcs(iv,imu,imu+1,iz,is) = +vfac*code_dt*vpa(iv-1)*mu(imu)*nux(iv-1,imu,iz)*mw(iv-1,imu,iz)/mw(iv-1,imu+1,iz)* dmu(imu-1)/dmu(imu) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
cc_blcs(iv,imu,imu-1,iz,is) = +vfac*code_dt*vpa(iv+1)*mu(imu)*nux(iv+1,imu,iz)*mw(iv+1,imu,iz)/mw(iv+1,imu-1,iz)* dmu(imu)/dmu(imu-1) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
cc_blcs(iv,imu,imu+1,iz,is) = -vfac*code_dt*vpa(iv+1)*mu(imu)*nux(iv+1,imu,iz)*mw(iv+1,imu,iz)/mw(iv+1,imu+1,iz)* dmu(imu-1)/dmu(imu) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
! mu operator
if (mu_operator) then
! quantities at mu_{i+1/2} and mu_{i-1/2}:
mup = 0.5*(mu(imu)+mu(imu+1))
mum = 0.5*(mu(imu)+mu(imu-1))
mwp = maxwell_vpa(iv,is)*exp(-2*bmag(ia,iz)*mup)
mwm = maxwell_vpa(iv,is)*exp(-2*bmag(ia,iz)*mum)
xp = sqrt(vpa(iv)**2 + 2*bmag(ia,iz)*mup)
xm = sqrt(vpa(iv)**2 + 2*bmag(ia,iz)*mum)
nuDp = spec(is)%vnew(is)*(erf(xp)-(erf(xp)-xp*(2/sqrt(pi))*exp(-xp**2)) / (2*xp**2)) / xp**3
nuDm = spec(is)%vnew(is)*(erf(xm)-(erf(xm)-xm*(2/sqrt(pi))*exp(-xm**2)) / (2*xm**2)) / xm**3
nupap = spec(is)%vnew(is)*2*(erf(xp)-xp*(2/sqrt(pi))*exp(-xp**2)) / (2*xp**2) / xp**3
nupam = spec(is)%vnew(is)*2*(erf(xm)-xm*(2/sqrt(pi))*exp(-xm**2)) / (2*xm**2) / xm**3
gam_mu = 2*(nupa(iv,imu,iz)*mu(imu)**2+nuD(iv,imu,iz)*vpa(iv)**2/(2*bmag(ia,iz))*mu(imu))*mw(iv,imu,iz)
gam_mum = 2*(nupam*mum**2+nuDm*vpa(iv)**2/(2*bmag(ia,iz))*mum)*mwm
gam_mup = 2*(nupap*mup**2+nuDp*vpa(iv)**2/(2*bmag(ia,iz))*mup)*mwp
! mu_operator (interior treatment):
bb_blcs(iv,imu,imu,ikxkyz) = bb_blcs(iv,imu,imu,ikxkyz) - code_dt*((gam_mu*(dmu(imu)/dmu(imu-1) - dmu(imu-1)/dmu(imu)) / (dmu(imu-1)+dmu(imu)) - gam_mum / dmu(imu-1))*dmu(imu)/dmu(imu-1) &
+(-gam_mu*(dmu(imu)/dmu(imu-1) - dmu(imu-1)/dmu(imu)) / (dmu(imu-1)+dmu(imu)) - gam_mup / dmu(imu))* dmu(imu-1)/dmu(imu)) / mw(iv,imu,iz) * 2./(dmu(imu-1)+dmu(imu))
bb_blcs(iv,imu,imu-1,ikxkyz) = bb_blcs(iv,imu,imu-1,ikxkyz) - code_dt*((-gam_mu * 1*dmu(imu)/dmu(imu-1) / (dmu(imu-1)+dmu(imu)) + gam_mum* 1/dmu(imu-1)) * dmu(imu)/dmu(imu-1) &
+gam_mu * 1*dmu(imu)/dmu(imu-1) / (dmu(imu-1)+dmu(imu)) * dmu(imu-1)/dmu(imu)) / mw(iv,imu-1,iz) * 2./(dmu(imu-1)+dmu(imu))
bb_blcs(iv,imu,imu+1,ikxkyz) = bb_blcs(iv,imu,imu+1,ikxkyz) - code_dt*( gam_mu * dmu(imu-1)/dmu(imu) / (dmu(imu-1)+dmu(imu)) * dmu(imu)/dmu(imu-1) &
+ (gam_mup/dmu(imu) - gam_mu*dmu(imu-1)/dmu(imu) / (dmu(imu-1)+dmu(imu))) * dmu(imu-1)/dmu(imu)) / mw(iv,imu+1,iz) * 2/(dmu(imu-1)+dmu(imu))
! mu operator, mixed (interior treatment):
aa_blcs(iv,imu,imu,iz,is) = aa_blcs(iv,imu,imu,iz,is) + code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv-1,imu,iz)*(dmu(imu)/dmu(imu-1)-dmu(imu-1)/dmu(imu))) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
aa_blcs(iv,imu,imu-1,iz,is) = aa_blcs(iv,imu,imu-1,iz,is) - code_dt*(vpa(iv)*mu(imu-1)*nux(iv,imu-1,iz)*mw(iv,imu-1,iz)*1/mw(iv-1,imu-1,iz)* dmu(imu)/dmu(imu-1)) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
aa_blcs(iv,imu,imu+1,iz,is) = aa_blcs(iv,imu,imu+1,iz,is) + code_dt*(vpa(iv)*mu(imu+1)*nux(iv,imu+1,iz)*mw(iv,imu+1,iz)*1/mw(iv-1,imu+1,iz)* dmu(imu-1)/dmu(imu)) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
cc_blcs(iv,imu,imu,iz,is) = cc_blcs(iv,imu,imu,iz,is) - code_dt*(vpa(iv)*mu(imu)*nux(iv,imu,iz)*mw(iv,imu,iz)*1/mw(iv+1,imu,iz)*(dmu(imu)/dmu(imu-1)-dmu(imu-1)/dmu(imu))) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
cc_blcs(iv,imu,imu-1,iz,is) = cc_blcs(iv,imu,imu-1,iz,is) + code_dt*(vpa(iv)*mu(imu-1)*nux(iv,imu-1,iz)*mw(iv,imu-1,iz)*1/mw(iv+1,imu-1,iz)* dmu(imu)/dmu(imu-1)) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
cc_blcs(iv,imu,imu+1,iz,is) = cc_blcs(iv,imu,imu+1,iz,is) - code_dt*(vpa(iv)*mu(imu+1)*nux(iv,imu+1,iz)*mw(iv,imu+1,iz)*1/mw(iv+1,imu+1,iz)* dmu(imu-1)/dmu(imu)) / (dmu(imu-1)+dmu(imu)) / (2*dvpa)
end if
end if
end if
! add gyro-diffusive term:
bb_blcs(iv,imu,imu,ikxkyz) = bb_blcs(iv,imu,imu,ikxkyz) + code_dt*cfac*0.5*kperp2(iky,ikx,ia,iz)*(spec(is)%smz/bmag(ia,iz))**2*(nupa(iv,imu,iz)*bmag(ia,iz)*mu(imu) + nuD(iv,imu,iz)*(vpa(iv)**2 + bmag(ia,iz)*mu(imu)))
end do
end do
end do
! switch from block storage to full matrix
blockmatrix = 0.
do ikxkyz = kxkyz_lo%llim_proc, kxkyz_lo%ulim_proc
iz = iz_idx(kxkyz_lo,ikxkyz)
is = is_idx(kxkyz_lo,ikxkyz)
do iv = 1, nvpa
! diagonal blocks:
blockmatrix(nmu*(iv-1)+1:nmu*iv, nmu*(iv-1)+1:nmu*iv, ikxkyz) = bb_blcs(iv,:,:,ikxkyz)
if (iv < nvpa) then
! subdiagonal blocks:
blockmatrix(nmu*iv+1:nmu*(iv+1), nmu*(iv-1)+1:nmu*iv, ikxkyz) = aa_blcs(iv+1,:,:,iz,is)
! superdiagonal blocks:
blockmatrix(nmu*(iv-1)+1:nmu*iv, nmu*iv+1:nmu*(iv+1), ikxkyz) = cc_blcs(iv,:,:,iz,is)
end if
end do
end do
! switch from full matrix to band-storage, for LAPACK banded solver routines:
! a_ij is stored in aband(ku+1+i-j,j) for $\max(1,j-ku) \leq i \leq \min(m,j+kl)$
cdiffmat_band = 0.
do ikxkyz = kxkyz_lo%llim_proc, kxkyz_lo%ulim_proc
iz = iz_idx(kxkyz_lo,ikxkyz)
is = is_idx(kxkyz_lo,ikxkyz)
do iv = 1, nmu*nvpa
do imu = 1, nmu*nvpa
if ((max(1, iv-(nmu+1)) .le. imu) .and. (imu .le. min(nvpa*nmu, iv+(nmu+1)))) then
cdiffmat_band(nmu+1 + nmu+1 + 1+imu-iv, iv, ikxkyz) = blockmatrix(imu, iv, ikxkyz)
end if
end do
end do
end do
! AVB: LU factorise cdiffmat, using LAPACK's zgbtrf routine for factorising banded matrices:
nc = nvpa*nmu
nb = nmu+1
lldab = 3*(nmu+1)+1
do ikxkyz = kxkyz_lo%llim_proc, kxkyz_lo%ulim_proc
call zgbtrf(nc, nc, nb, nb, cdiffmat_band(:,:,ikxkyz), lldab, ipiv(:,ikxkyz), info)
end do
end subroutine init_fp_diffmatrix
subroutine init_fp_conserve
use finite_differences, only: tridag
use linear_solve, only: lu_decomposition
use stella_time, only: code_dt
use species, only: nspec, spec
use zgrid, only: nzgrid, ntubes
use vpamu_grids, only: ztmax, maxwell_vpa, maxwell_mu
use vpamu_grids, only: nmu, nvpa,vpa, vperp2
use vpamu_grids, only: set_vpa_weights, wgts_vpa, wgts_mu
use kt_grids, only: naky, nakx
use stella_layouts, only: kxkyz_lo
use stella_layouts, only: iky_idx, ikx_idx, iz_idx, is_idx, it_idx
use dist_fn_arrays, only: gvmu
use gyro_averages, only: aj0v
use fields, only: get_fields, get_fields_by_spec
use stella_geometry, only: bmag
use job_manage, only: time_message, timer_local
implicit none
integer :: ikxkyz, iky, ikx, iz, is, it
integer :: imu, idx, ix, iv
logical :: conservative_wgts
real :: dum2, dum3
complex, dimension (:,:,:,:), allocatable :: dum1
complex, dimension (:,:,:,:,:), allocatable :: field
complex, dimension (:,:), allocatable :: sumdelta
complex, dimension (:,:), allocatable :: gvmutr
complex, dimension (:), allocatable :: ghrs
if (.not.allocated(vpadiff_response)) then
nresponse_vpa = 1
if (momentum_conservation) nresponse_vpa = nresponse_vpa + nspec
if (energy_conservation) nresponse_vpa = nresponse_vpa + nspec
allocate (vpadiff_response(nresponse_vpa,nresponse_vpa,kxkyz_lo%llim_proc:kxkyz_lo%ulim_alloc))
vpadiff_response = 0.
allocate (vpadiff_idx(nresponse_vpa,kxkyz_lo%llim_proc:kxkyz_lo%ulim_alloc))
end if
allocate (dum1(naky,nakx,-nzgrid:nzgrid,ntubes))
allocate (field(naky,nakx,-nzgrid:nzgrid,ntubes,nspec))
allocate (sumdelta(nmu,nvpa)); sumdelta = 0.
allocate (gvmutr(nvpa,nmu))
allocate (ghrs(nmu*nvpa))
! set wgts to be equally spaced to ensure exact conservation properties
conservative_wgts = .true.
call set_vpa_weights (conservative_wgts)
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)
is = is_idx(kxkyz_lo,ikxkyz)
! AVB: calculate Green's function: supply unit impulse to rhs at location (iv,imu), solve for response
sumdelta = 0.
do iv = 1, nvpa
do imu = 1, nmu
ghrs = 0.
ghrs(nmu*(iv-1)+imu) = 1.
! solve for response:
call zgbtrs('No transpose', nvpa*nmu, nmu+1, nmu+1, 1, cdiffmat_band(:,:,ikxkyz), 3*(nmu+1)+1, ipiv(:,ikxkyz), ghrs, nvpa*nmu, info)
do ix = 1, nvpa
gvmutr(ix,:) = ghrs(nmu*(ix-1)+1 : nmu*(ix-1)+nmu)
end do
sumdelta = sumdelta + ztmax(iv,is)*maxwell_mu(1,iz,imu,is)*aj0v(imu,ikxkyz)*transpose(gvmutr)
end do
end do
gvmu(:,:,ikxkyz) = transpose(sumdelta) ! AVB: gvmu(:,:,ikxkyz) now contains: sum_iv sum_imu [ response(vpar,mu)_{iv,imu} ] / phi^(n+1) = h_{hom}^n+1(vpa,mu)_ikxkyz / phi^{n+1}_ikxkyz
end do
! gvmu contains dhs/dphi
! for phi equation, need 1-P[dhs/dphi]
call get_fields (gvmu, field(:,:,:,:,1), dum1, dist='h')
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)
vpadiff_response(1,1,ikxkyz) = 1.0-field(iky,ikx,iz,it,1) ! AVB: ravelled response
end do
! now get LU decomposition for vpadiff_response
do ikxkyz = kxkyz_lo%llim_proc, kxkyz_lo%ulim_proc
call lu_decomposition (vpadiff_response(:,:,ikxkyz),vpadiff_idx(:,ikxkyz),dum2)
end do
! reset wgts to default setting
conservative_wgts = .false.
call set_vpa_weights (conservative_wgts)
deallocate (dum1, field)
end subroutine init_fp_conserve
subroutine init_vpadiff_matrix
use stella_time, only: code_dt
use species, only: nspec, spec
use vpamu_grids, only: dvpa, vpa, nvpa
use stella_layouts, only: kxkyz_lo
use stella_layouts, only: iky_idx, ikx_idx, iz_idx, is_idx
use stella_geometry, only: bmag
use dist_fn_arrays, only: kperp2
implicit none
integer :: ikxkyz, iky, ikx, iz, is
integer :: ia
if (.not.allocated(aa_vpa)) allocate (aa_vpa(nvpa,nspec))
! if (.not.allocated(bb_vpa)) allocate (bb_vpa(nvpa,nspec))
if (.not.allocated(bb_vpa)) allocate (bb_vpa(nvpa,kxkyz_lo%llim_proc:kxkyz_lo%ulim_alloc))
if (.not.allocated(cc_vpa)) allocate (cc_vpa(nvpa,nspec))
! deal with boundary points (BC is f(vpa)=0 beyond +/- vpa_max)
aa_vpa(1,:) = 0.0 ; cc_vpa(nvpa,:) = 0.0
! 2nd order centered differences for d/dvpa (1/2 dh/dvpa + vpa h)
do is = 1, nspec
aa_vpa(2:,is) = -code_dt*spec(is)%vnew(is)*0.5*(1.0/dvpa-vpa(:nvpa-1))/dvpa
! bb_vpa(:,is) = 1.0+code_dt*spec(is)%vnew(is)/dvpa**2
cc_vpa(:nvpa-1,is) = -code_dt*spec(is)%vnew(is)*0.5*(1.0/dvpa+vpa(2:))/dvpa
end do
ia = 1
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)
is = is_idx(kxkyz_lo,ikxkyz)
bb_vpa(:,ikxkyz) = 1.0 + code_dt*spec(is)%vnew(is) &
* (0.25*kperp2(iky,ikx,ia,iz)*(spec(is)%smz/bmag(ia,iz))**2 + 1./dvpa**2)
end do
end subroutine init_vpadiff_matrix
subroutine init_mudiff_matrix
use stella_time, only: code_dt
use species, only: nspec, spec
use zgrid, only: nzgrid
use stella_geometry, only: bmag
use vpamu_grids, only: dmu, mu, nmu
use stella_layouts, only: kxkyz_lo
use stella_layouts, only: iky_idx, ikx_idx, iz_idx, is_idx
use dist_fn_arrays, only: kperp2
implicit none
integer :: ikxkyz, iky, ikx, iz, is
integer :: ia
! TMP FOR TESTING -- MAB
! integer :: imu
real, dimension (:), allocatable :: dmu_ghost, dmu_cell, mu_cell
! add ghost cell at mu=0 and beyond mu_max for purposes of differentiation
! note assuming here that grid spacing for ghost cell is equal to
! grid spacing for last non-ghost cell
allocate (dmu_ghost(nmu))
dmu_ghost(:nmu-1) = dmu ; dmu_ghost(nmu) = dmu(nmu-1)
! this is mu at cell centres (including to left and right of mu grid boundary points)
allocate (mu_cell(nmu))
mu_cell(:nmu-1) = 0.5*(mu(:nmu-1)+mu(2:))
mu_cell(nmu) = mu(nmu)+0.5*dmu(nmu-1)
! this is mu_{j+1/2} - mu_{j-1/2}
allocate (dmu_cell(nmu))
dmu_cell(1) = mu_cell(1)
dmu_cell(2:) = mu_cell(2:)-mu_cell(:nmu-1)
if (.not.allocated(aa_mu)) allocate (aa_mu(-nzgrid:nzgrid,nmu,nspec))
if (.not.allocated(bb_mu)) allocate (bb_mu(nmu,kxkyz_lo%llim_proc:kxkyz_lo%ulim_alloc))
if (.not.allocated(cc_mu)) allocate (cc_mu(-nzgrid:nzgrid,nmu,nspec))
ia = 1
! deal with boundary points (BC is f(mu)=0 beyond mu_max and collision operator vanishes for mu -> 0)
aa_mu(:,1,:) = 0.0 ; cc_mu(:,nmu,:) = 0.0
! 2nd order centered differences for dt * nu * d/dmu (mu/B*dh/dmu + 2*mu*h)
do is = 1, nspec
do iz = -nzgrid, nzgrid
aa_mu(iz,2:,is) = -code_dt*spec(is)%vnew(is)*mu_cell(:nmu-1)*(1.0/(bmag(ia,iz)*dmu)-1.0)/dmu_cell(2:)
cc_mu(iz,:nmu-1,is) = -code_dt*spec(is)%vnew(is)*mu_cell(:nmu-1)*(1.0/(bmag(ia,iz)*dmu)+1.0)/dmu_cell(:nmu-1)
end do
end do
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)
is = is_idx(kxkyz_lo,ikxkyz)
bb_mu(1,ikxkyz) = 1.0 + code_dt*spec(is)%vnew(is) &
* (0.25*kperp2(iky,ikx,ia,iz)*(spec(is)%smz/bmag(ia,iz))**2 &
+ 1.0/(dmu(1)*bmag(ia,iz)) - 1.0)
bb_mu(2:nmu-1,ikxkyz) = 1.0 + code_dt*spec(is)%vnew(is) &
* (0.25*kperp2(iky,ikx,ia,iz)*(spec(is)%smz/bmag(ia,iz))**2 &
+ (mu_cell(2:nmu-1)/dmu(2:)+mu_cell(:nmu-2)/dmu(:nmu-2)) &
/(dmu_cell(2:nmu-1)*bmag(ia,iz)) - 1.0)
bb_mu(nmu,ikxkyz) = 1.0 + code_dt*spec(is)%vnew(is) &
* (0.25*kperp2(iky,ikx,ia,iz)*(spec(is)%smz/bmag(ia,iz))**2 &
+ mu_cell(nmu-1)*(1.0/(dmu(nmu-1)*bmag(ia,iz)) + 1.0)/dmu_cell(nmu))
end do
deallocate (dmu_ghost, dmu_cell, mu_cell)
end subroutine init_mudiff_matrix
subroutine init_vpadiff_conserve
use finite_differences, only: tridag
use linear_solve, only: lu_decomposition
use stella_time, only: code_dt
use species, only: nspec, spec
use zgrid, only: nzgrid, ntubes
use vpamu_grids, only: ztmax, maxwell_vpa, maxwell_mu
use vpamu_grids, only: nmu, vpa, vperp2
use vpamu_grids, only: set_vpa_weights
use kt_grids, only: naky, nakx
use stella_layouts, only: kxkyz_lo
use stella_layouts, only: iky_idx, ikx_idx, iz_idx, it_idx, is_idx
use dist_fn_arrays, only: gvmu
use gyro_averages, only: aj0v
use fields, only: get_fields, get_fields_by_spec
implicit none
integer :: ikxkyz, iky, ikx, iz, it, is
integer :: imu
integer :: idx
logical :: conservative_wgts
real :: dum2
complex, dimension (:,:,:,:), allocatable :: dum1
complex, dimension (:,:,:,:,:), allocatable :: field
if (.not.allocated(vpadiff_response)) then
nresponse_vpa = 1
if (momentum_conservation) nresponse_vpa = nresponse_vpa + nspec
if (energy_conservation) nresponse_vpa = nresponse_vpa + nspec
allocate (vpadiff_response(nresponse_vpa,nresponse_vpa,kxkyz_lo%llim_proc:kxkyz_lo%ulim_alloc))
vpadiff_response = 0.
allocate (vpadiff_idx(nresponse_vpa,kxkyz_lo%llim_proc:kxkyz_lo%ulim_alloc))
end if
allocate (dum1(naky,nakx,-nzgrid:nzgrid,ntubes))
allocate (field(naky,nakx,-nzgrid:nzgrid,ntubes,nspec))
! set wgts to be equally spaced to ensure exact conservation properties
conservative_wgts = .true.
call set_vpa_weights (conservative_wgts)
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)
is = is_idx(kxkyz_lo,ikxkyz)
do imu = 1, nmu
gvmu(:,imu,ikxkyz) = ztmax(:,is)*maxwell_mu(1,iz,imu,is)*aj0v(imu,ikxkyz)
call tridag (1, aa_vpa(:,is), bb_vpa(:,ikxkyz), cc_vpa(:,is), gvmu(:,imu,ikxkyz))
end do
end do
! gvmu contains dhs/dphi
! for phi equation, need 1-P[dhs/dphi]
! for upar equations, need -Us[dhs/dphi]
! for energy conservation, need -Qs[dhs/dphi]
!!! FLAG DSO - The following lines are NOT appropriate for
!!! zonal modes with adianbatic electrons!
call get_fields (gvmu, field(:,:,:,:,1), dum1, dist='h')
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)
vpadiff_response(1,1,ikxkyz) = 1.0-field(iky,ikx,iz,it,1)
end do
idx = 2
if (momentum_conservation) then
call get_upar (gvmu, field)
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)
vpadiff_response(idx:idx+nspec-1,1,ikxkyz) = -field(iky,ikx,iz,it,:)
end do
idx = idx + nspec
end if
if (energy_conservation) then
call get_temp (gvmu, field)
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)
vpadiff_response(idx:idx+nspec-1,1,ikxkyz) = -field(iky,ikx,iz,it,:)
end do
end if
idx = 2
if (momentum_conservation) then
do ikxkyz = kxkyz_lo%llim_proc, kxkyz_lo%ulim_proc
iz = iz_idx(kxkyz_lo,ikxkyz)
is = is_idx(kxkyz_lo,ikxkyz)
do imu = 1, nmu