-
Notifications
You must be signed in to change notification settings - Fork 3
/
m_multislice.f90
914 lines (718 loc) · 35.6 KB
/
m_multislice.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
!--------------------------------------------------------------------------------
!
! Copyright (C) 2017 L. J. Allen, H. G. Brown, A. J. D’Alfonso, S.D. Findlay, B. D. Forbes
!
! This program is free software: you can redistribute it and/or modify
! it under the terms of the GNU General Public License as published by
! the Free Software Foundation, either version 3 of the License, or
! (at your option) any later version.
!
! This program is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU General Public License for more details.
!
! You should have received a copy of the GNU General Public License
! along with this program. If not, see <http://www.gnu.org/licenses/>.
!
!--------------------------------------------------------------------------------
module m_multislice
use m_precision, only: fp_kind
implicit none
logical :: save_grates = .false.
logical :: load_grates = .false.
character(1024) :: grates_filename
integer*4::qep_mode
logical :: output_probe_intensity = .false.,additional_transmission_function = .false.,pure_phase
logical,allocatable :: output_cell_list(:)
real(fp_kind),allocatable :: output_thickness_list(:)
integer,allocatable :: cell_map(:)
character*200,allocatable::amplitude_fnam(:),phase_fnam(:)
complex(fp_kind),allocatable,dimension(:,:) :: shift_arrayx, shift_arrayy
integer(4) :: n_slices !number of subslices in the unit cell
integer(4) :: maxnat_slice !maximum number of atoms for a particular type in the supercell
real(fp_kind), allocatable :: a0_slice(:,:)
integer(4), allocatable :: nat_slice(:,:) !number of atoms for each type in the slice (supercell)
integer(4), allocatable :: nat_slice_unitcell(:,:) !number of atoms for each type in the slice (unit cell)
real(fp_kind), allocatable :: tau_slice(:,:,:,:) !the atom positions for each potential subslice (supercell)
real(fp_kind), allocatable :: tau_slice_unitcell(:,:,:,:) !the atom positions for each potential subslice (unit cell)
real(fp_kind), allocatable :: prop_distance(:) !Unit cell subslice propagation distance
real(fp_kind), allocatable :: depths(:) !Unit cell subslice depths (i.e. same as above)
real(fp_kind), allocatable :: ss_slice(:,:)
!DMH
#ifdef LIN
character(len=1), parameter :: path_sep='/'
#elif WIN
character(len=1), parameter :: path_sep='\'
#else
#error "path_sep not defined, Set this constant for your system"
#endif
interface load_save_add_grates
module procedure load_save_add_grates_qep,load_save_add_grates_abs
end interface
contains
!This subroutine samples from the available phase grates and then performs one iteration of the multislice algorithm (called in CPU versions only)
subroutine qep_multislice_iteration(psi,propagator,transmission,nopiy,nopix,ifactory,ifactorx,idum,n_qep_grates,mode,shift_arrayy,shift_arrayx)
use m_numerical_tools,only:ran1
complex(fp_kind),intent(inout)::psi(nopiy,nopix)
complex(fp_kind),intent(in)::propagator(nopiy,nopix),transmission(nopiy,nopix,n_qep_grates),shift_arrayy(nopiy),shift_arrayx(nopix)
integer*4,intent(in)::nopiy,nopix,n_qep_grates,mode,ifactory,ifactorx
integer*4,intent(inout)::idum
integer*4::shifty,shiftx,nran
complex(fp_kind)::trans(nopiy,nopix)
! Phase grate
nran = floor(n_qep_grates*ran1(idum)) + 1
if(mode == 1) then !On the fly calculation
!call make_qep_potential(trans, tau_slice, nat_slice, ss_slice(7,j))
!psi_out = psi*trans
elseif(mode == 2) then !Quick shift
shiftx = floor(ifactorx*ran1(idum)) * nopix/ifactorx
shifty = floor(ifactory*ran1(idum)) * nopiy/ifactory
trans = cshift(cshift(transmission(:,:,nran),shifty,dim=1),shiftx,dim=2)
call multislice_iteration(psi,propagator,trans,nopiy,nopix)
elseif(mode == 3) then !Phase ramp shift
shiftx = floor(ifactorx*ran1(idum)) + 1
shifty = floor(ifactory*ran1(idum)) + 1
call phase_shift_array(transmission(:,:,nran),trans,shift_arrayy,shift_arrayx)
call multislice_iteration(psi,propagator,trans,nopiy,nopix)
else
call multislice_iteration(psi,propagator,transmission(:,:,nran),nopiy,nopix)
endif
end subroutine
!This subroutine performs one iteration of the multislice algorithm (called in CPU versions only)
!Probe (psi) input and output is in real space
subroutine multislice_iteration(psi,propagator,transmission,nopiy,nopix)
use cufft_wrapper
complex(fp_kind),intent(inout)::psi(nopiy,nopix)
complex(fp_kind),intent(in)::propagator(nopiy,nopix),transmission(nopiy,nopix)
integer*4,intent(in)::nopiy,nopix
! Transmit through slice potential
psi = psi*transmission
! Propagate to next slice
call fft2(nopiy,nopix,psi,nopiy,psi,nopiy)
psi = psi*propagator
call ifft2(nopiy,nopix,psi,nopiy,psi,nopiy)
end subroutine
subroutine prompt_output_probe_intensity
use m_user_input, only: get_input
use global_variables, only: thickness, nopiy, nopix
use output, only: output_prefix,split_filepath
use m_string, only: to_string,command_line_title_box
implicit none
integer :: i_output,j
real(fp_kind) :: thickness_interval
character(1024)::dir,fnam
logical::invalid_thickness
call command_line_title_box('Output probe intensity')
write(*,*) 'The probe intensity as a function of thickness'
write(*,*) 'can be outputted to file at each probe position.'
write(*,*) 'The user is advised that the outputted dataset'
write(*,*) 'may be very large.'
write(*,*)
10 write(*,*) '<0> Proceed'
write(*,*) '<1> Output probe intensity'
call get_input('<0> Proceed <1> Output probe intensity', i_output)
write(*,*)
output_probe_intensity = i_output ==1
if(output_probe_intensity) then
invalid_thickness = .true.
do while(invalid_thickness)
write(*,*) 'At what thickness interval (in Angstroms)',char(10),' should intensities be outputted?'
call get_input('At what thickness interval should intensities be outputted?', thickness_interval)
write(*,*)
invalid_thickness = thickness_interval.le.0.0_fp_kind .or. thickness_interval.gt.thickness
if(invalid_thickness) write(*,*) 'ERROR: invalid thickness.'
enddo
call split_filepath(output_prefix,dir,fnam)
!DMH call system('mkdir '//trim(adjustl(dir))//'\Probe_intensity')
call system('mkdir '//trim(adjustl(dir))//path_sep//'Probe_intensity')
call generate_cell_list(thickness_interval)
call write_thicknesss_to_file
write(*,*) 'The probe intensities will be written to the files'
write(*,*)
!DMH write(*,*) ' '//trim(adjustl(dir))//'\Probe_intensity' // trim(adjustl(fnam)) // '_ProbeIntensity*.bin'
write(*,*) ' '//trim(adjustl(dir))//path_sep//'Probe_intensity' // trim(adjustl(fnam)) // '_ProbeIntensity*.bin'
write(*,*)
if (fp_kind.eq.4) write(*,*) 'as 32-bit big-endian floats.'
if (fp_kind.eq.8) write(*,*) 'as 64-bit big-endian floats.'
write(*,*) 'Each file contains a sequence of ' // to_string(nopiy) // 'x' // to_string(nopix) // ' arrays.'
write(*,*)
end if
end subroutine
subroutine generate_cell_list(thickness_interval)
use global_variables, only: thickness, ncells, a0
implicit none
real(fp_kind) :: thickness_interval
integer :: i_cell, count,i,j
real(fp_kind) :: t,tout
if(allocated(output_cell_list)) deallocate(output_cell_list)
allocate(output_cell_list(maxval(ncells)*n_slices))
output_cell_list = .false.
if(allocated(cell_map)) deallocate(cell_map)
allocate(cell_map(maxval(ncells)*n_slices))
cell_map = 0
t = 0.0_fp_kind
tout = thickness_interval
count = 0
do i= 1,maxval(ncells)
do j=1,n_slices
if((t+depths(j+1)*a0(3).gt.tout)) then
output_cell_list((i-1)*n_slices+j) =.true.
count = count + 1
tout = (floor((t+depths(j+1)*a0(3))/thickness_interval)+1)*thickness_interval
endif
enddo
t = t+a0(3)
enddo
if(allocated(output_thickness_list)) deallocate(output_thickness_list)
allocate(output_thickness_list(count))
t = 0.0_fp_kind
tout = thickness_interval
count = 0
do i= 1,maxval(ncells)
do j=1,n_slices
if((t+depths(j+1)*a0(3).gt.tout)) then
count = count + 1
output_thickness_list(count) = t+depths(j+1)*a0(3)
cell_map((i-1)*n_slices+j) = count
tout = (floor((t+depths(j+1)*a0(3))/thickness_interval)+1)*thickness_interval
endif
enddo
t = t+a0(3)
enddo
end subroutine
subroutine write_thicknesss_to_file
use output, only: output_prefix
implicit none
integer :: i,j
character(1024) :: filename,dir,fnam
j = index(output_prefix,'/',back=.true.)
j = max(j,index(output_prefix,'\',back=.true.))
if(j>0) then
dir = trim(adjustl(output_prefix(:j)))
fnam = trim(adjustl(output_prefix(j:)))
!DMH filename = trim(adjustl(dir))//'\Probe_intensity'//trim(adjustl(fnam))//'_probe_intensity_thicknesss.txt'
filename = trim(adjustl(dir))//path_sep//'Probe_intensity'//trim(adjustl(fnam))//'_probe_intensity_thicknesss.txt'
else
!DMH filename = 'Probe_intensity\'//trim(adjustl(output_prefix))//'_probe_intensity_thicknesss.txt'
filename = 'Probe_intensity'//path_sep//trim(adjustl(output_prefix))//'_probe_intensity_thicknesss.txt'
endif
write(*,*) 'The thicknesses at which the probe intensity'
write(*,*) 'is being outputted have been written to'
write(*,*)
write(*,*) ' ' // trim(filename)
write(*,*)
open(unit=8734, file=filename)
do i = 1, size(output_thickness_list)
write(8734, *) output_thickness_list(i)
enddo
close(8734)
end subroutine
subroutine probe_intensity_to_file(probe_intensity,i_df,ny,nx,n_qep_passes,probe_ndf,nysample,nxsample)
use output, only: output_prefix,quad_shift
use global_variables, only: nopiy,nopix
use m_string
real(fp_kind),intent(in)::probe_intensity(nopiy,nopix,size(output_thickness_list))
integer*4,intent(in)::i_df,ny,nx,n_qep_passes,probe_ndf,nysample,nxsample
integer*4::j,z
character*1024::filename,fnam,dir
j = index(output_prefix,'/',back=.true.)
j = max(j,index(output_prefix,'\',back=.true.))
z = size(output_thickness_list)
if(j>0) then
dir = trim(adjustl(output_prefix(:j)))
fnam = trim(adjustl(output_prefix(j:)))
!DMH filename = trim(adjustl(dir))//'\Probe_intensity\'//trim(adjustl(fnam))//'_ProbeIntensity'
filename = trim(adjustl(dir))//path_sep//'Probe_intensity\'//trim(adjustl(fnam))//'_ProbeIntensity'
else
!DMH filename = 'Probe_intensity\'//trim(adjustl(output_prefix))//'_ProbeIntensity'
filename = 'Probe_intensity'//path_sep//trim(adjustl(output_prefix))//'_ProbeIntensity'
endif
if (probe_ndf.gt.1) filename = trim(filename) // '_df' // to_string(i_df)
if (nysample.gt.1) filename = trim(filename) // '_ny' // to_string(ny)
if (nxsample.gt.1) filename = trim(filename) // '_nx' // to_string(nx)
filename = trim(filename) // '_'//to_string(nopiy)//'x'//to_string(nopix)//'x'//to_string(z)//'.bin'
open(4985, file=filename, form='binary', convert='big_endian')
do j=1,z
write(4985) quad_shift(probe_intensity(:,:,j),nopiy,nopix)/ n_qep_passes
enddo
close(4985)
end subroutine
subroutine prompt_save_load_grates
use m_user_input, only: get_input
use output, only: output_prefix
use m_string
implicit none
integer :: i_save_load, i_retry,i
logical :: exists,retry
i_save_load = -1
do while(i_save_load<0.or.i_save_load>2)
call command_line_title_box('Save/load transmission functions')
write(*,*) 'Warning: the files outputted when saving may be very large.'
write(*,*)
write(*,*) '<0> Proceed without saving or loading'
write(*,*) '<1> Save transmission functions'
write(*,*) '<2> Load transmission functions'
write(*,*) '<3> Add additional transmission function '
write(*,*) ' (eg. from magnetic structure) from file'
call get_input('<0> continue <1> save <2> load', i_save_load)
write(*,*)
323 format( ' - The xtl file',/,' - The slicing of the unit cell',/,&
&' - The choice of thermal scattering model (QEP vs. absorptive)',/,&
&' - The tiling of the unit cell',/,' - The number of pixels',/,&
&' - (For absorptive model: whether absorption is included)',/,&
&' - (For QEP model: the number of distinct transmission functions)',/,&
&' - (For QEP model: phase ramp shift choice)',/)
select case (i_save_load)
case (0)
return
case (1)
save_grates = .true.
grates_filename = trim(adjustl(output_prefix)) // '_transmission_functions.bin'
write(*,*) 'The transmission functions will be saved to the file'
write(*,*)
write(*,*) ' ' // trim(grates_filename)
write(*,*)
write(*,*) 'They can be loaded for later calculations provided'
write(*,*) 'the following parameters are identical:'
write(6,323)
case (2)
write(*,*) 'It is up to the user to ensure that the parameters used'
write(*,*) 'to create the loaded transmission functions are consistent'
write(*,*) 'with those of the current calculation:'
write(6,323)
retry=.true.
do while(retry)
write(*,*) 'Enter filename of transmission functions:'
call get_input('filename of transmission functions', grates_filename)
write(*,*)
inquire(file=grates_filename, exist=exists)
retry=.not.exists
if (.not.exists) then
write(*,*) 'ERROR: cannot find this file.'
write(*,*) '<1> Enter again'
write(*,*) '<2> Proceed without loading'
call get_input('<1> Enter again <2> Proceed without loading', i_retry)
write(*,*)
retry = i_retry==1
!DMH if(.not.retry) return
if(.not.retry) then
return
endif
endif
enddo
load_grates = .true.
case (3)
additional_transmission_function=.true.
pure_phase=.true.
!DMH if(.not.allocated(amplitude_fnam)) allocate(amplitude_fnam(n_slices),phase_fnam(n_slices))
if(.not.allocated(amplitude_fnam)) then
allocate(amplitude_fnam(n_slices),phase_fnam(n_slices))
endif
do i=1,n_slices
if(.not.pure_phase) then
write(*,*) char(10),' Please input filename of amplitude of additional transmission function for slice ',i
if(i==1) then
write(*,*) 'If your intended transmission function is a pure phase object (ie. the '
write(*,*) 'amplitude of theadditional transmission function is everywhere 1), input'
write(*,*) '-1 and no amplitude file willbe loaded nor will you be prompted for the '
write(*,*) 'amplitude for the remaining slices'
endif
call get_input('Amplitude of additional transmission function',amplitude_fnam(i))
pure_phase = trim(adjustl(amplitude_fnam(i)))=='-1'
endif
write(*,*) char(10),' Please input filename of phase of additional transmission function for slice ',i
call get_input('Phase of additional transmission function',phase_fnam(i))
enddo
end select
enddo
end subroutine
subroutine load_save_add_grates_qep(idum,qep_grates,nopiy,nopix,n_qep_grates,n_slices,nt,nat_slice)
use m_numerical_tools, only: gasdev
use output
use global_variables, only: ak,pi
integer(4),intent(inout):: idum
complex(fp_kind),intent(inout)::qep_grates(nopiy,nopix,n_qep_grates,n_slices)
integer*4,intent(in)::nopiy,nopix,n_qep_grates,n_slices,nt,nat_slice(nt,n_slices)
integer*4::j,i,m,n
real(fp_kind)::junk
real(fp_kind)::amplitude(nopiy,nopix),phase(nopiy,nopix)
if (load_grates) then
write(*,*) 'Loading transmission functions from file...'
write(*,*)
open(unit=3984, file=grates_filename, form='binary', convert='big_endian')
read(3984) qep_grates
close(3984)
! Make sure the random number sequence is as if grates were calculated
! So call gasdev as many times as it would usually be called
do j = 1, n_slices;do i = 1, n_qep_grates;do m = 1, nt;do n = 1, nat_slice(m,j)*2
junk = gasdev(idum)
enddo;enddo;enddo;enddo
return
endif
if(additional_transmission_function) then
write(*,*) 'Adding addition transmission function to file...'
amplitude = 1
do j=1,n_slices
!DMH if(.not.pure_phase) call binary_in(nopiy,nopix,amplitude,amplitude_fnam(j))
if(.not.pure_phase) then
call binary_in(nopiy,nopix,amplitude,amplitude_fnam(j))
endif
call binary_in(nopiy,nopix,phase,phase_fnam(j))
qep_grates(:,:,:,j) = qep_grates(:,:,:,j)+spread(transpose(phase),dim=3,ncopies=n_qep_grates)/pi/a0_slice(3,j)*ak
if(.not.pure_phase) then
call binary_in(nopiy,nopix,amplitude,amplitude_fnam(j))
qep_grates(:,:,:,j) = qep_grates(:,:,:,j)+spread(cmplx(0_fp_kind,log(transpose(amplitude))),dim=3,ncopies=n_qep_grates)
endif
enddo
qep_mode=4
endif
if (save_grates) then
write(*,*) 'Saving transmission functions to file...',char(10)
open(unit=3984, file=grates_filename, form='binary', convert='big_endian')
write(3984) qep_grates*pi*spread(spread(spread(a0_slice(3,1:n_slices),dim=1,ncopies=nopiy),dim=2,ncopies=nopix),dim=3,ncopies=n_qep_grates)/ak;close(3984)
endif
end subroutine
subroutine load_save_add_grates_abs(abs_grates,nopiy,nopix,n_slices)
use global_variables, only: pi,ak
use m_numerical_tools, only: gasdev
use output
complex(fp_kind),intent(inout)::abs_grates(nopiy,nopix,n_slices)
integer*4,intent(in)::nopiy,nopix,n_slices
integer*4::j,i,m,n
real(fp_kind)::junk
real(fp_kind)::amplitude(nopiy,nopix),phase(nopiy,nopix)
if (load_grates) then
write(*,*) 'Loading transmission functions from file...'
write(*,*)
open(unit=3984, file=grates_filename, form='binary', convert='big_endian')
read(3984) abs_grates
close(3984)
return
endif
if(additional_transmission_function) then
write(*,*) 'Adding addition transmission function from file...'
amplitude = 1
do j=1,n_slices
call binary_in(nopiy,nopix,phase,phase_fnam(j))
abs_grates(:,:,j) = abs_grates(:,:,j)+transpose(phase)/pi/a0_slice(3,j)*ak
if(.not.pure_phase) then
call binary_in(nopiy,nopix,amplitude,amplitude_fnam(j))
abs_grates(:,:,j) = abs_grates(:,:,j)+cmplx(0_fp_kind,log(transpose(amplitude)))
endif
enddo
endif
if (save_grates) then
write(*,*) 'Saving transmission functions to file...',char(10)
open(unit=3984, file=grates_filename, form='binary', convert='big_endian')
write(3984) abs_grates*pi*spread(spread(a0_slice(3,1:n_slices),dim=1,ncopies=nopiy),dim=2,ncopies=nopix)/ak;close(3984)
endif
end subroutine
function make_detector(nopiy,nopix,ifactory,ifactorx,ss,kmin,kmax,phi,delphi) result(detector)
use m_crystallography
!makes a detector on an array nopiy x nopix with Fourier space pixel size deltaky x deltakx
!which measures from kmin to kmax, supplying phi and delphi (detector orientation
!and angular range in radians) will create a detector segment
integer*4,intent(in)::nopiy,nopix,ifactory,ifactorx
real(fp_kind),intent(in)::kmin,kmax,phi,delphi,ss(7)
optional::phi,delphi
real(fp_kind)::phi_(2),g_vec_array(3,nopiy,nopix),ang,kabs,detector(nopiy,nopix),pi
integer*4::y,x
logical::segment,notwrapped
pi = 4.0d0*atan(1.0d0)
!If the phi and delphi variables are present, then make a segmented detector
segment = present(phi).and.present(delphi)
if(segment) then
!Force angle to be between 0 and 2pi
phi_= mod([phi,phi+delphi],2*pi)
do y=1,2
if (phi_(y)<0) phi_(y) = phi_(y) + 2*pi
enddo
notwrapped = phi_(2)>phi_(1)
endif
detector = 0
call make_g_vec_array(g_vec_array,ifactory,ifactorx)
detector = 0
do y=1,nopiy;do x=1,nopix
kabs = trimr(g_vec_array(:,y,x),ss)
if((kabs.ge.kmin).and.(kabs.le.kmax)) then
detector(y,x) =1
if(segment) then
!Force angle to be between 0 and 2pi
ang = modulo(atan2(g_vec_array(1,y,x),g_vec_array(2,y,x)),2*pi)
if (notwrapped) then
!DMH if(.not.(ang>phi_(1).and.ang<=phi_(2))) detector(y,x)=0
if(.not.(ang>phi_(1).and.ang<=phi_(2))) then
detector(y,x)=0
endif
else
!DMH if(.not.(ang>phi_(1).or.ang<=phi_(2))) detector(y,x)=0
if(.not.(ang>phi_(1).or.ang<=phi_(2))) then
detector(y,x)=0
endif
endif
endif;endif
enddo;enddo;
end function
subroutine setup_qep_parameters(n_qep_grates,n_qep_passes,nran,quick_shift,ifactory,ifactorx)
use m_user_input, only: get_input
use m_string, only: to_string,command_line_title_box
implicit none
integer*4,intent(out)::n_qep_grates,n_qep_passes,nran
integer*4,intent(in)::ifactory,ifactorx
logical,intent(in)::quick_shift
integer :: i
call command_line_title_box('QEP model parameters')
qep_mode=1
do while(qep_mode<2.or.qep_mode>4)
10 write(*,*) 'Enter the number of distinct transmission functions to calculate:'
call get_input("Number of phase grates calculated", n_qep_grates)
write(*,*)
write(*,*) 'As implemented in muSTEM, the QEP multislice algorithm will randomly shift the'
write(*,*) 'unit cells in your supercell around to effectively generate new random '
write(*,*) 'transmission functions. This means the QEP calculation samples a larger '
write(*,*) 'effective number of random transmission functions than would otherwise be the'
write(*,*) 'case.'
write(*,*)
write(*,*) 'If your choice of unit cell tiling and pixel grid size is such that each unit '
write(*,*) 'cell is an integer number of pixels, unit cell shifting can be achieved by '
write(*,*) 'circular shift of the transmission function arrays in memory (quick shifting).'
write(*,*) 'Otherwise the unit cells will have to be shifted with sub-pixel precision '
write(*,*) 'using the Fourier shift algorithm which will impact calculation speed. '
write(*,*)
write(*,*) 'A third option is to not shift the unit cells around at all, but the user is'
write(*,*) 'advised that this requires a much larger number of distinct transmission '
write(*,*) 'functions to achieve convergence than would usually be the case.'
write(*,*)
if (quick_shift) then
write(6,100) to_string(ifactorx*ifactory), to_string(n_qep_grates), to_string(ifactorx*ifactory*n_qep_grates)
100 format(' The choice of tiling and grid size permits quick shifting.',/,&
&' The effective number of transmission functions used in ',/,&
&' calculations will be ', a, ' * ', a, ' = ', a, '.', /)
qep_mode=2
else
write(6,101)
101 format( ' Your choice of tiling and grid size does not permit quick shifting ', /, &
&' of the precalculated transmission functions. Shifting using the ', /, &
&' Fourier shift algorithm can be performed but is time consuming. ', /, &
&' You may wish to go back and calculate more distinct transmission ', /, &
&' functions, or simply proceed without using phase ramp shifting. ' /)
i=-1
do while(i<1.or.i>3)
110 write(6,111)
111 format( ' <1> Go back and choose a larger number', /, &
&' <2> Proceed with phase ramp shifting', /, &
&' <3> Proceed without phase ramp shifting', / &
)
call get_input("<1> choose more <2> phase ramp <3> no phase ramp", i)
write(*,*)
if (i.eq.2 .and. (ifactory.gt.1 .or. ifactorx.gt.1)) then
qep_mode=3
call setup_phase_ramp_shifts
elseif (i.eq.3) then
qep_mode=4
endif
enddo
endif
enddo
write(6,*) 'Enter the number of passes to perform for QEP calculation:'
write(*,*) 'Warning: using only a single pass is usually NOT sufficient.'
call get_input("Number of Monte Carlo calculated", n_qep_passes )
write(*,*)
write(6,*) 'Enter the starting position of the random number sequence:'
call get_input("Number of ran1 discarded", nran )
write(*,*)
end subroutine
subroutine setup_phase_ramp_shifts
use global_variables, only: nopiy, nopix, ifactory, ifactorx
implicit none
integer(4) :: i
real(fp_kind) :: r_coord
if(allocated(shift_arrayy)) deallocate(shift_arrayy)
if(allocated(shift_arrayx)) deallocate(shift_arrayx)
allocate(shift_arrayy(nopiy,ifactory))
allocate(shift_arrayx(nopix,ifactorx))
do i = 1,ifactory
r_coord=float(i)/float(ifactory)
call make_shift_oned(shift_arrayy(:,i),nopiy,r_coord)
enddo
do i = 1,ifactorx
r_coord=float(i)/float(ifactorx)
call make_shift_oned(shift_arrayx(:,i),nopix,r_coord)
enddo
end subroutine
subroutine setup_slicing_depths()
use m_user_input, only: get_input
use m_precision, only: fp_kind
use m_string,only:command_line_title_box
implicit none
integer :: i_slice
call command_line_title_box('Unit cell slicing')
22 write(6,23) char(143)
23 format(' Do you wish to slice the unit cell in the beam direction?', /, &
&' This may be a good idea if the unit cell is larger than 2 ', a1, /, &
&' in the z-direction.', /, &
&' <1> Yes <2> No ')
call get_input('Slice unit cell <1> yes <2> no', i_slice)
write(*,*)
if (i_slice.eq.1) then
call calculate_depths_slicing
elseif (i_slice.eq.2) then
call calculate_depths_no_slicing
else
goto 22
endif
end subroutine
subroutine calculate_depths_no_slicing
implicit none
n_slices = 1
allocate(depths(2))
depths(1) = 0.0_fp_kind
depths(2) = 1.0_fp_kind
end subroutine
subroutine calculate_depths_slicing
use m_user_input
implicit none
integer(4) i, ichoice
n_slices = -1
do while(n_slices<1)
write(*,*) 'Enter the number of slices per unit cell in the beam direction:'
call get_input("Number of distinct potentials ", n_slices)
write(*,*)
if (n_slices.eq.1) then; call calculate_depths_no_slicing
elseif (n_slices.gt.1) then
if (allocated(depths)) deallocate(depths)
allocate(depths(n_slices+1))
depths(n_slices+1) = 1.0_fp_kind
write(6,10)
10 format( ' You will now be asked to enter the depths at which slicing',/,&
&' will take place. These should be entered as fractions of',/,&
&' the unit cell. It is the front of the slice which should',/,&
&' be entered. (e.g. three even slices would be entered as',/,&
&' 0.0, 0.3333, 0.6666 ).', /)
ichoice = 0
do while(ichoice<1.or.ichoice>2)
write(6,16)
16 format(' How would you like to specify the slice depths?', /, &
&' <1> Manually ', /, &
&' <2> Automatically (uniformly spaced)')
call get_input("<1> manual or <2> auto slicing", ichoice)
write(*,*)
if(ichoice.eq.1) then
! Manual slicing
do i = 1, n_slices
write(6,20,ADVANCE='NO') achar(13), i
!DMH
call flush(6)
20 format( a1,' Enter the fractional depth of slice number ', i4, ': ')
call get_input("depths", depths(i))
enddo
elseif (ichoice.eq.2) then
! Automatic slicing
21 format( ' Fractional depth of slice ', i4, ': ', f7.4)
do i = 1, n_slices
depths(i) = float( (i-1)) / float(n_slices)
write(6,21) i, depths(i)
enddo
endif
enddo
endif
enddo
write(*,*)
end subroutine
subroutine calculate_slices
use global_variables, only: nat, ifactory, ifactorx, nt, a0, deg, tau, nm
use m_crystallography, only: cryst
implicit none
integer :: j
maxnat_slice = maxval(nat) * ifactory * ifactorx
allocate(nat_slice(nt,n_slices),tau_slice(3,nt,maxnat_slice,n_slices),&
&a0_slice(3,n_slices),ss_slice(7,n_slices),prop_distance(n_slices),&
nat_slice_unitcell(nt,n_slices),tau_slice_unitcell(3,nt,nm,n_slices))
do j = 1, n_slices
a0_slice(1,j) = a0(1)*ifactorx
a0_slice(2,j) = a0(2)*ifactory
if (j .eq. n_slices) then
a0_slice(3,j) = (1.0_fp_kind-depths(j))*a0(3)
else
a0_slice(3,j) = (depths(j+1)-depths(j))*a0(3)
endif
call cryst(a0_slice(:,j),deg,ss_slice(:,j))
! Make supercell cell subsliced
call calculate_tau_nat_slice( depths(j),depths(j+1),tau,tau_slice(:,:,:,j),nm,nt,nat,nat_slice(:,j),ifactorx,ifactory )
! Make unit cell subsliced
call make_mod_tau_unitcell( depths(j),depths(j+1),tau,tau_slice_unitcell(:,:,:,j),nm,nt,nat,nat_slice_unitcell(:,j) )
enddo
prop_distance = a0_slice(3,:)
end subroutine
subroutine calculate_tau_nat_slice( depth1, depth2, tau, tau_slice, nm, nt, nat, nat_slice, ifactorx, ifactory )
implicit none
integer(4) nm, nt, ifactorx, ifactory
integer(4) nat(nt)
integer(4) nat_slice(nt)
real(fp_kind) depth2, depth1, tau(3,nt,nm)
real(fp_kind) tau_slice(3,nt,nm*ifactorx*ifactory)
integer(4) i,j,jj, mm, nn
real(fp_kind) diff, tol
real(fp_kind) factorx, factory
tol = 1.0d-4
diff = depth2 - depth1
factorx = float( ifactorx )
factory = float( ifactory )
do i = 1, nt
jj = 1
do mm = 1, ifactory
do nn = 1, ifactorx
do j = 1, nat(i)
if( (tau(3,i,j) .lt. (depth2-tol)) .and.(tau(3,i,j) .ge. (depth1-tol)) ) then
tau_slice(1,i,jj) = tau(1,i,j)/factorx + float(nn-1)/factorx
tau_slice(2,i,jj) = tau(2,i,j)/factory + float(mm-1)/factory
tau_slice(3,i,jj) = (tau(3,i,j)-depth1)/diff
jj = jj + 1
elseif(diff.eq.1.0) then
tau_slice(1,i,jj) = tau(1,i,j)/factorx + float(nn-1)/factorx
tau_slice(2,i,jj) = tau(2,i,j)/factory + float(mm-1)/factory
tau_slice(3,i,jj) = (tau(3,i,j)-depth1)/diff
endif
enddo
enddo
enddo
nat_slice(i) = jj-1
enddo
end subroutine
!--------------------------------------------------------------------------------------
subroutine make_mod_tau_unitcell( depth1, depth2, tau, mod_tau, nm, nt, nat, nat2)
!subroutine make_mod_tau_force() is called if natural depths
!were NOT chosen. In this case the current and next depth are
!required to select which of the tau fit within this slice.
!Note that a tolerance is employed -- this measures how far
!beyond the limits an atom can be and still be counted.
!CAREFUL to realise that false holz can be an issue using this forced
!slicing.
use m_precision
implicit none
integer(4) nm, nt
integer(4) nat(nt),nat2(nt)
real(fp_kind) depth2, depth1, tau(3,nt,nm),mod_tau(3,nt,nm)
integer(4) i,j,jj
real(fp_kind) diff, tol
tol = 1.0e-4_fp_kind
diff = depth2-depth1
do i=1,nt
jj=1
do j=1,nat(i)
if( (tau(3,i,j) .lt. (depth2-tol)) .and.(tau(3,i,j) .ge. (depth1-tol)) ) then
mod_tau(1,i,jj) = tau(1,i,j)
mod_tau(2,i,jj) = tau(2,i,j)
mod_tau(3,i,jj) = (tau(3,i,j)-depth1)/diff
jj = jj+1
elseif(diff.eq.1.0) then
mod_tau(1,i,jj) = tau(1,i,j)
mod_tau(2,i,jj) = tau(2,i,j)
mod_tau(3,i,jj) = (tau(3,i,j)-depth1)/diff
endif
enddo
nat2(i)=jj-1
enddo
end subroutine
end module