-
Notifications
You must be signed in to change notification settings - Fork 2
/
crystal3.f90
executable file
·3488 lines (3029 loc) · 119 KB
/
crystal3.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
!
!
! TODO:
! - check sqrt(psi_h*psi_hbar) in the formulas, in many cases is cdabs(psi_h)
! - check asymmetry factor in ML and PP (not appearing?)
! - extinction depth is in intensity whereas pendellosung depth is in amplitude
! - check the sign of c in ML Bragg
!
!
!----
!----
!---- MODULE: crystal3
!----
!---- Main module for crystal diffraction
!----
!----
!----
!----
!----
!---- The sequence of calling is the following
!----
!----
!---- crystal_loadcrystaldata: to load preprocessor file (now written by XOP)
!---- crystal_set: to set (align) the crystal at a given energy
!---- crystal_info:
!---- crysta_{prefect,mosaic,bent_pp,bent_ml} to compute reflectivity
!----
!----
!---- the interface is the single "crystal" routine that should be called with
!---- "kwhat" flag:
!---- <0 (calls crystal_loadcrystaldata)
!---- 0 (calls crystal_set)
!---- 1 (calls crystal_info)
!---- 2 gets reflectivity
!----
!----
Module crystal3
!---- Use Modules ----!
use shadow_globaldefinitions
use shadow_math
use elasticity
!---- Variables ----!
implicit none
!public
!todo move to a common place consistenly with shadow
real(kind=8) :: r0 = 2.817940289458D-13 ! classical electron radius in cm
!todo remove
!xianbos b_factor in c calculation (12.96 for ml model)(-0.0077 for pp model)
real(kind=8) :: B_factor
!
! type (structure) to store crystallographic data from file
! (data usually created by a preprocessor)
!
integer(kind=ski),public,parameter :: NMAXENER=1000, NMAXATOMS=100
!type, public, bind(C) :: crystalData
type, public :: crystalData
!C RN: the constant (e^2/mc^2)/V or the ration between the classical e-
!C radius and the volume of the unit cell [cm^(-2)]
real(kind=skr) :: rn,d_spacing
integer(kind=ski) :: itype
real(kind=skr),dimension(NMAXATOMS) :: atnum
real(kind=skr),dimension(NMAXATOMS) :: temper
integer(kind=ski) :: nbatom
complex(kind=skx),dimension(NMAXATOMS) :: G,G_BAR
real(kind=skr),dimension(NMAXATOMS) :: G_0
integer(kind=ski) :: nf0coeff
real(kind=skr),dimension(NMAXATOMS,11) :: f0coeff
integer(kind=ski) :: npoint
real(kind=skr),dimension(NMAXENER) :: energy
real(kind=skr),dimension(NMAXATOMS,NMAXENER) :: fp,fpp,fcompton
real(kind=skr),dimension(NMAXATOMS) :: fract
character(len=sklen) :: filename
real(kind=skr) :: version
end type crystalData
!
! this is the type that contains all the crystal information needed
! by crystal calculations.
!
!type, public, bind(C) :: oeSetup
type, public :: oeSetup
! these are user inputs
real(kind=skr) :: thickness=0.1D0 ! crystal thickness
integer(kind=ski) :: f_bragg_a=0 ! asymmetric?
real(kind=skr) :: a_bragg=0.0D0 ! asymmetry angle
integer(kind=ski) :: f_mosaic=0 ! [0] perfect, [1]mosaic,
! [2] bent multilamellar
! [3] bent Penning-Polder
real(kind=skr) :: spread_mos=0.0D0 ! crystal mosaicity
integer(kind=ski) :: f_refrac=0 ! geometry: 0=Braggref, 1=Laueref
! 2=Braggtrans, 3=Lauetrans
character(len=sklen) :: file_refl ! file created by preprocessor
! this stores the data in the preprocessor file
type(crystalData) :: crystalData ! type to load xtal info from file_refl
! this is defined in module: elasticity, stores the inputs and
! calculated values
type(crystalElasticity) :: crystalElasticity ! type to load elastic info
! this are inputs specific to bent crystals
Real(kind=skr) :: bent_rm=0d0 !meridional bending radius Rm
Real(kind=skr) :: bent_rs=0d0 !Sagital bending radius Rs [cm]
! these are calculated variables (from crystal_set). Bragg angles
Real(kind=skr) :: GRAZE=0d0 !Bragg angle
Real(kind=skr) :: THETA_B=0d0 !corrected incidence angle,
Real(kind=skr) :: THETA_B_H=0d0 !corrected reflected angle
Real(kind=skr) :: THETA_B_SYM=0d0 !corrected Bragg angle for sym crystal
Real(kind=skr) :: THETA_B_COR=0d0 !corrected Bragg angle
Real(kind=skr) :: L_EXT_S=0d0 ! Pendellosung length s-pol
Real(kind=skr) :: L_EXT_P=0d0 ! idem p-pol
Real(kind=skr) :: ASS_FAC=0d0 ! asymmetry factor
Real(kind=skr) :: SSR=0d0 ! half-Darwin width s-pol
Real(kind=skr) :: SPR=0d0 ! idem p-plo
Real(kind=skr) :: Q_MOS=0d0 ! mosaics Q value
! bent crystal calculations (e.g., bent_rs will contain the anticlastic)
Real(kind=skr) :: bent_m1=0d0 !meridional bending moment
Real(kind=skr) :: bent_m2=0d0 !sagittal bending moment
Real(kind=skr) :: bent_rm_calc=0d0 !meridional bending radius Rm
Real(kind=skr) :: bent_rs_calc=0d0 !Sagital bending radius Rs [cm]
Real(kind=skr) :: ml_c=0d0 ! c(sigma) parameter of the ML
!Real(kind=skr) :: pp_beta=0d0 ! beta parameters iof the PP
!Real(kind=skr) :: pp_betac=0d0 ! beta-critical parameters iof the PP
Real(kind=skr) :: pp_G=0d0 ! G parameter of the PP
! these are the 2pi/lambda and photon energy at which the crystal has
! been set (by crystal_set)
Real(kind=skr) :: Q_PHOT_SET=0d0 ,PHOT_SET=0d0
!vectors
real(kind=skr),dimension(3) :: VNOR ! normal to the surface (pointing out)
real(kind=skr),dimension(3) :: VPAR ! on the crystal, in diffraction plane
real(kind=skr),dimension(3) :: AXIS ! on the crystal, perp to diff plane
real(kind=skr),dimension(3) :: BH ! crystal reciprocal lattice vector
real(kind=skr),dimension(3) :: BH0 ! normalized crystal reciprocal lattice vector
real(kind=skr),dimension(3) :: VIN_BRAGG ! Bragg direction
! (corrected for refraction)
real(kind=skr),dimension(3) :: VIN_BRAGG_UNCORR ! Bragg direction
! (not corrected for refraction)
real(kind=skr),dimension(3) :: VIN_BRAGG_ENERGY
! and their respective output vectors
real(kind=skr),dimension(3) :: VOUT_BRAGG
real(kind=skr),dimension(3) :: VOUT_BRAGG_UNCORR
real(kind=skr),dimension(3) :: VOUT_BRAGG_ENERGY
end type oeSetup
!---- Everything is private unless explicitly made public ----!
private
!---- List of public routines ----!
public :: crystal, scat_direction
public :: crystal_loadCrystalData, plotgle, plotidl, crystal_set_vectors
!---- List of private routines ----!
private :: crystal_printCrystalData
private :: CRYSTAL_FH, CRYSTAL_INFO, CRYSTAL_SET
private :: CRYSTAL_PERFECT, CRYSTAL_MOSAIC,CRYSTAL_BENT_PP,CRYSTAL_BENT_ML
private :: CRYSTAL_BENT_MOMENTS,CRYSTAL_BENT_PP_CALC_BETA,CRYSTAL_BENT_ML_CALC_C
Contains
!
!---- Routines ----!
!
!C
!C subroutine scat_direction
!C calculates the output scattering direction applying the
!C equations k_out_par = k_in_par + BH_par and |k_in| = |k_out|
!C
!C kk% = ( 2 pi / lambda ) * vv%
!C vv% direction versor
!C q_phot = 2 pi / lambda [cm**-1]
!C bh = ( 2 pi / d_spacing ) normal_to_bragg_planes
!C vnor = surface_normal
!C
subroutine scat_direction (vvin1,bh,vnor1,q_phot,vvout)
implicit none
real(kind=skr),dimension(3),intent(in) :: vvin1,bh,vnor1
real(kind=skr),intent(in) :: q_phot
real(kind=skr),dimension(3),intent(out) :: vvout
real(kind=skr) :: mod2_kkout_par,mod_kkout_perp,tmp1,tmp2
real(kind=skr),dimension(3) :: vnor,vvin,kkin,kkin_par,kkin_perp
real(kind=skr),dimension(3) :: kkout,kkout_par,kkout_perp
real(kind=skr),dimension(3) :: bh_par,bh_perp
real(kind=skr),dimension(3) :: vvout1,vvout2,vvoutEwald
vnor = vnor1
vvin = vvin1
Call norm (vnor,vnor)
Call norm (vvin,vvin)
Call scalar (vvin,q_phot,kkin)
Call proj (kkin,vnor,kkin_perp)
Call vector (kkin_perp,kkin,kkin_par)
Call proj (bh,vnor,bh_perp)
Call vector (bh_perp,bh,bh_par)
Call vsum (kkin_par,bh_par,kkout_par)
Call dot(kkout_par,kkout_par,mod2_kkout_par )
mod_kkout_perp = dsqrt(q_phot**2 - mod2_kkout_par)
Call scalar (vnor,mod_kkout_perp,kkout_perp)
! plus sign in mod_kkout_perp
Call vsum (kkout_par,kkout_perp,kkout)
Call norm (kkout,vvout1)
! minus sign in mod_kkout_perp
Call vsum (kkout_par,-kkout_perp,kkout)
Call norm (kkout,vvout2)
! select vvout1 or vvout2
vvoutEwald = bh+vvin*q_phot
call norm(vvoutEwald,vvoutEwald)
call dot(vvout1,vvoutEwald,tmp1)
call dot(vvout2,vvoutEwald,tmp2)
tmp1 = 1.0 - abs(tmp1)
tmp2 = 1.0 - abs(tmp2)
vvout = vvout1
if (tmp2.lt.tmp1) vvout = vvout2
return
end subroutine scat_direction
!C+++
!C SUBROUTINE CRYSTAL_FH
!C
!C PURPOSE Computes the structute factor of a crystal
!C from data stored in crystalData type
!C
!C ALGORITHM
!C
!C MODIFIED Created by M. Sanchez del Rio (Feb. 1996)
!C
!C Modified by MSR 96/12/19 to include the number
!C of coefficients in the f0 parametrization. This
!C number is 9 for the CromerMann data and 11 for
!C the Waasmaier Kirfel data. Note that the crystal
!C file has changed!
!C
!C---
!C
!C INPUT PARAMETERS:
!C fUnit: fileUnit
!C fUnit = [-Inf,3) Not used
!C fUnit >= 3 Writes info in file unit=fUnit
!C oe1: oeSetup type with crystal data
!C PHOT [real]: Photon energy in eV
!C THETA [real]: Scattering grazing angle in rads
!C
!C OUTPUT PARAMETERS:
!C FH |
!C FH_BAR |
!C F_0 | Complex with the returned structure factors
!C PSI_H | F, Psi and the refraction index.
!C PSI_HBAR |
!C PSI_0 |
!C REFRAC |
!C
!C
SUBROUTINE CRYSTAL_FH(fUnit, oe1, PHOT, THETA, & ! inputs
FH,FH_BAR,F_0,PSI_H,PSI_HBAR,PSI_0,REFRAC,STRUCT) ! outputs
implicit none
integer(kind=ski),intent(in):: fUnit
type(oeSetup),intent(in) :: oe1
real(kind=skr),intent(in) :: PHOT,THETA
complex(kind=skx),intent(out) :: FH,FH_BAR,F_0
complex(kind=skx),intent(out) :: psi_h,psi_hbar,psi_0,REFRAC,STRUCT
integer(kind=ski):: if0center
integer(kind=ski):: i,j,nener,ierr
real(kind=skr),dimension(NMAXATOMS) :: F0,F1,F2,F2compton
real(kind=skr) :: r_lam0,ratio,version
real(kind=skr) :: delta_ref,absorp,absorpNoCompton
complex(kind=skx):: psi_conj,F_0compton
complex(kind=skx):: ci,F(NMAXATOMS)
integer(kind=ski):: i_debug=0
CI = (0.0D0,1.0D0)
if (i_debug.GT.0) then
write(i_debug,*) '<><> CRYSTAL_FH: fUnit=',fUnit
write(i_debug,*) '<><> CRYSTAL_FH: oe1%crystalData%d_spacing=',oe1%crystalData%d_spacing
write(i_debug,*) '<><> CRYSTAL_FH: oe1%crystalData%rn=',oe1%crystalData%rn
endif
!C
!C Computes structure factor at given wavelength and angle.
!C
if (i_debug.EQ.1) then
write (i_debug,*) '<><>CRYSTAL_FH: working at energy: ',phot
write (i_debug,*) '<><>CRYSTAL_FH: working at angle [rads] : ',theta
write (i_debug,*) '<><>CRYSTAL_FH: working at angle [deg] : ',theta*todeg
endif
IF (PHOT.LT.oe1%crystalData%ENERGY(1).OR.PHOT.GT.oe1%crystalData%ENERGY(oe1%crystalData%npoint)) THEN
write (*,*) 'CRYSTAL_FH: Error: Incoming photon energy is out of range: ',PHOT
RETURN
END IF
!C
!C Build the fo scattering form factor from its coefficients
!C
ratio = sin(theta)/(toangs/phot)
if (ratio.GT.2) then
write (*,*) 'CRYSTAL_FH: ratio sin(theta)/lambda > 2'
write (*,*) ' ratio : ',ratio
write (*,*) ' Paramatrizatiog for Fo may fail.'
end if
if (i_debug.gt.0) write (i_debug,*) 'CRYSTAL_FH: ratio is : ',ratio
if0center = (1+oe1%crystalData%nf0coeff)/2
do j=1,oe1%crystalData%nbatom
f0 (j) = oe1%crystalData%f0coeff(j,if0center)
do i=1,if0center-1
f0(j) = f0(j) + oe1%crystalData%f0coeff(j,i) * &
dexp(-1.0d0*oe1%crystalData%f0coeff(j,i+if0center)*ratio**2)
if (i_debug.gt.0) write (i_debug,*) '<><>CRYSTAL_FH: f0 i,j = ',f0(j),i,j
end do
if (i_debug.gt.0) write (i_debug,*) '<><>CRYSTAL_FH: F0(j) = ',F0(j)
end do
!C
!C Interpolate for the atomic scattering factor.
!C
DO 299 I = 1, oe1%crystalData%npoint
299 IF (oe1%crystalData%ENERGY(I).GT.PHOT) GO TO 101
!C
I = oe1%crystalData%npoint
101 NENER = I - 1
do j=1,oe1%crystalData%nbatom
F1(j)= oe1%crystalData%FP(j,NENER) + (oe1%crystalData%FP(j,NENER+1) - &
oe1%crystalData%FP(j,NENER)) * (PHOT - oe1%crystalData%ENERGY(NENER)) / &
(oe1%crystalData%ENERGY(NENER+1) - oe1%crystalData%ENERGY(NENER))
F2(j)= oe1%crystalData%FPP(j,NENER) + (oe1%crystalData%FPP(j,NENER+1) - &
oe1%crystalData%FPP(j,NENER)) * (PHOT - oe1%crystalData%ENERGY(NENER)) / &
(oe1%crystalData%ENERGY(NENER+1) - oe1%crystalData%ENERGY(NENER))
F2compton(j)= &
oe1%crystalData%FPP(j,NENER)*oe1%crystalData%Fcompton(j,NENER) +&
(oe1%crystalData%FPP(j,NENER+1)*oe1%crystalData%Fcompton(j,NENER+1) - &
oe1%crystalData%FPP(j,NENER)*oe1%crystalData%Fcompton(j,NENER)) * &
(PHOT - oe1%crystalData%ENERGY(NENER)) / &
(oe1%crystalData%ENERGY(NENER+1) - oe1%crystalData%ENERGY(NENER))
if (i_debug.gt.0) then
write (i_debug,*) '<><>CRYSTAL_FH: F1 = ',F1(j)
write (i_debug,*) '<><>CRYSTAL_FH: F2 = ',F2(j)
write (i_debug,*) '<><>CRYSTAL_FH: Fcompton(j,NENER) = ', &
oe1%crystalData%Fcompton(j,NENER)
endif
end do
R_LAM0 = TOCM/PHOT
if (i_debug.gt.0) write (i_debug,*) '<><>CRYSTAL_FH: R_LAM0 =', R_LAM0
do j=1,oe1%crystalData%nbatom
F(j)= F0(j) + F1(j) + CI*F2(j)
if (i_debug.gt.0) write (i_debug,*) '<><>CRYSTAL_FH: F = ',F(j)
end do
!C
!C FH and FH_BAR are the structure factors for (h,k,l) and (-h,-k,-l).
!C
F_0 = (0.0D0, 0.0D0)
F_0compton = (0.0D0, 0.0D0)
FH = (0.0D0, 0.0D0)
FH_BAR = (0.0D0, 0.0D0)
do i=1,oe1%crystalData%nbatom
FH = FH + (oe1%crystalData%G(i) * F(i) * oe1%crystalData%TEMPER(i) * &
oe1%crystalData%FRACT(i))
!C From Brennan's private communication:
!C Changing the temperature doesn't change F0, but it does change
!C Chi0, due to changing the size of the unit cell, which changes
!C gamma.
!C MSR 96/02/14
!
! todo: check whether the compton correction affects also f1.
! Now, it only affects the imaginary part f2 but not for the structure
! factor, only for the refraction index (and attenuation coeff).
F_0 = F_0 + oe1%crystalData%G_0(i) * (oe1%crystalData%atnum(i) + F1(i) + &
CI * F2(i)) * oe1%crystalData%FRACT(i)
F_0compton = F_0compton + &
oe1%crystalData%G_0(i) * (oe1%crystalData%atnum(i) + F1(i) + &
CI * F2compton(i)) * oe1%crystalData%FRACT(i)
FH_BAR = FH_BAR + (oe1%crystalData%G_BAR(i) * F(i) * &
oe1%crystalData%TEMPER(i) * oe1%crystalData%FRACT(i))
end do
STRUCT = sqrt(FH * FH_BAR)
if (i_debug.gt.0) then
write (i_debug,*) '<><>CRYSTAL_FH: FH = ',FH
write (i_debug,*) '<><>CRYSTAL_FH: FH_BAR = ',FH_BAR
write (i_debug,*) '<><>CRYSTAL_FH: f_0 = ',f_0
write (i_debug,*) '<><>CRYSTAL_FH: STRUCT = ',STRUCT
write (i_debug,*) ' '
write (i_debug,*) '<><>CRYSTAL_FH: oe1%crystalData%RN = ',oe1%crystalData%RN
write (i_debug,*) '<><>CRYSTAL_FH: r_lam0 = ',r_lam0
endif
!C
!C PSI_CONJ = F*( note: PSI_HBAR is PSI at -H position and is
!C proportional to fh_bar but PSI_CONJ is complex conjugate os PSI_H)
!C
! srio@esrf.eu 20140116 added sign (see Zachariasen pag 114, eq. 3.95)
psi_h = -oe1%crystalData%rn*r_lam0**2/pi*fh
psi_hbar = -oe1%crystalData%rn*r_lam0**2/pi*fh_bar
psi_0 = -oe1%crystalData%rn*r_lam0**2/pi*f_0
psi_conj = -oe1%crystalData%rn*r_lam0**2/pi*conjg(fh)
if (i_debug.gt.0) then
write (i_debug,*) '<><>CRYSTAL_FH: PSI_H = ',PSI_H
write (i_debug,*) '<><>CRYSTAL_FH: PSI_HBAR = ',PSI_HBAR
write (i_debug,*) '<><>CRYSTAL_FH: PSI_0 = ',PSI_0
endif
!C
!C computes refractive index.
!C ([3.171] of Zachariasen's book)
!C
!
! values without compton correction
!
REFRAC = (1.0D0,0.0D0) - R_LAM0**2*oe1%crystalData%RN*F_0/TWOPI
DELTA_REF = 1.0D0 - DREAL(REFRAC)
ABSORP = 2.0d0 * TWOPI *(-DIMAG(REFRAC)) / R_LAM0
if (i_debug.gt.0) then
write (i_debug,*) '<><>CRYSTAL_FH: REFRAC NO_COMPTON= ', REFRAC
write (i_debug,*) '<><>CRYSTAL_FH: DELTA_REF NO_COMPTON= ', DELTA_REF
write (i_debug,*) '<><>CRYSTAL_FH: ABSORP NO_COMPTON= ', ABSORP
endif
!
! values WITH compton correction
!
absorpNoCompton = absorp
REFRAC = (1.0D0,0.0D0) - R_LAM0**2*oe1%crystalData%RN*F_0compton/TWOPI
DELTA_REF = 1.0D0 - DREAL(REFRAC)
ABSORP = 2.0d0 * TWOPI *(-DIMAG(REFRAC)) / R_LAM0
if (i_debug.EQ.1) then
write (i_debug,*) '<><>CRYSTAL_FH: REFRAC = ', REFRAC
write (i_debug,*) '<><>CRYSTAL_FH: DELTA_REF = ', DELTA_REF
write (i_debug,*) '<><>CRYSTAL_FH: ABSORP = ', ABSORP
endif
!C
!C if fUnit > 3, write text to unit fUnit
!C
if (fUnit.gt.3) then
write(fUnit,*) ' '
write(fUnit,*) '******************************************************'
write(fUnit,*) ' ** DIFF_PAT v1.8 (24 Mar 2014) ** '
write(fUnit,*) 'Crystal data from file '//trim(oe1%crystalData%filename)
WRITE(fUnit,*) ' at energy = ',PHOT,' eV'
WRITE(fUnit,*) ' = ',R_LAM0*1E8,' Angstroms'
WRITE(fUnit,*) ' and at angle = ',THETA*TODEG,' degrees'
WRITE(fUnit,*) ' = ',THETA,' rads'
write(fUnit,*) '******************************************************'
WRITE(fUnit,*) ' '
do j=1,oe1%crystalData%nbatom
WRITE(fUnit,*) 'For atom ',j,':'
WRITE(fUnit,*) ' fo + fp+ i fpp = '
WRITE(fUnit,*) ' ',f0(j),'+',f1(j),'+ i',f2(j),'= '
WRITE(fUnit,*) ' ',f0(j)+f1(j)+ci*f2(j)
WRITE(fUnit,*) ' Z = ',oe1%crystalData%atnum(j)
end do
WRITE(fUnit,*) 'Structure factor F(0,0,0) = ',F_0
WRITE(fUnit,*) 'Structure factor FH = ',FH
WRITE(fUnit,*) 'Structure factor FH_BAR = ',FH_BAR
WRITE(fUnit,*) 'Structure factor F(h,k,l) = sqrt(FH*FH_BAR) ',STRUCT
WRITE(fUnit,*) 'Psi_0 = ',psi_0
WRITE(fUnit,*) 'Psi_H = ',psi_h
WRITE(fUnit,*) 'Psi_HBar = ',psi_hbar
WRITE(fUnit,*) '|Psi_H/Psi_HBar| = ',cdabs(psi_h/psi_hbar)
WRITE(fUnit,*) 'sqrt(Psi_H*Psi_HBar) = ',sqrt(psi_h*psi_hbar)
WRITE(fUnit,*) '|Psi_H| = ',cdabs(psi_h)
WRITE(fUnit,*) 'Refraction index = 1 - delta - i*beta :'
WRITE(fUnit,*) ' delta = ',DELTA_REF
WRITE(fUnit,*) ' beta = ',-1.0D0*DIMAG(REFRAC)
WRITE(fUnit,*) 'Absorption coeff = ',ABSORP,' cm^-1'
WRITE(fUnit,*) 'Absorption correction (mu_total/mu_photoelec) = ',ABSORP/absorpNoCompton
WRITE(fUnit,*) ' '
WRITE(fUnit,*) 'e^2/(mc^2)/V = ',oe1%crystalData%rn,' cm^-2'
WRITE(fUnit,*) 'd-spacing = ',oe1%crystalData%d_spacing*1.0E8,' Angstroms'
WRITE(fUnit,*) 'SIN(theta)/Lambda = ',Ratio
WRITE(fUnit,*) ' '
endif
END SUBROUTINE CRYSTAL_FH
!C+++
!C SUBROUTINE crystal_loadCrystalData
!C
!C PURPOSE Stores the crystallographic data from externam file
!C created with a preprocessor (xop/xcrystal or bragg)
!C
!C ALGORITHM
!C
!C MODIFIED Created by M. Sanchez del Rio
!C
!C Modified by MSR 96/12/19 to include the number
!C of coefficients in the f0 parametrization. This
!C number is 9 for the CromerMann data and 11 for
!C the Waasmaier Kirfel data. Note that the crystal
!C file has changed!
!C Modified by MSR 2011/08/16 to port to g95 and
!C integrate in Shadow3
!C
!C---
!C
!C INPUT AND OUTPUT PARAMETERS:
!C FILEIN : input, the name of the file with data
!C xtal: output, the type containing the data loaded from file.
!C
!C
SUBROUTINE crystal_loadCrystalData (FILEIN,xtal)
implicit none
character(len=sklen),intent(in) :: FILEIN
type(crystalData),intent(out) :: xtal
integer(kind=ski) :: if0center
integer(kind=ski) :: i,j,nener,ierr
integer(kind=ski) :: i_debug=0
character(len=sklen):: text
if (i_debug.gt.0) write(i_debug,*) '<><> crystal_loadCrystalData: file is: '//trim(FILEIN)
OPEN (25,FILE=FILEIN,STATUS='OLD', FORM='FORMATTED',ERR=77)
read(25,'(A)',err=79) text
read(25,*,err=79) xtal%version, xtal%itype
if (xtal%itype.NE.1) then
write (*,*) 'CRYSTAL_LOADCRYSTALDATA: Error: Data file type not yet implemented.'
stop
end if
read(25,'(A)',err=79) text
read(25,*,err=79) xtal%RN,xtal%d_spacing !,TEMPER
read(25,'(A)',err=79) text
read(25,*,err=79) xtal%nbatom
if (xtal%nbatom.GT.NMAXATOMS) then
write (*,*) 'CRYSTAL_LOADCRYSTALDATA: Error: Maximum number of atoms allowad: ',NMAXATOMS
end if
read(25,'(A)',err=79) text
read(25,*,err=79) (xtal%ATNUM(I),i=1,xtal%nbatom)
read(25,'(A)',err=79) text
read(25,*,err=79) (xtal%FRACT(I),i=1,xtal%nbatom)
read(25,'(A)',err=79) text
read(25,*,err=79) (xtal%TEMPER(I),i=1,xtal%nbatom)
read(25,'(A)',err=79) text
read(25,*,err=79) (xtal%G_0(i),i=1,xtal%nbatom)
read(25,'(A)',err=79) text
do i=1,xtal%nbatom
read(25,*,err=79) xtal%G(I)
read(25,*,err=79) xtal%G_BAR(I)
end do
read(25,'(A)',err=79) text
do i=1,xtal%nbatom
read(25,*,err=79) xtal%nf0coeff,(xtal%f0coeff(i,j),j=1,xtal%nf0coeff)
end do
read(25,'(A)',err=79) text
read(25,*,err=79) xtal%NPOINT
read(25,'(A)',err=79) text
do I = 1, xtal%NPOINT
read(25,*,err=79) xtal%energy(i)
if ( (xtal%version-2.39).GT.0) then
! new version (since 2.4) implements compton correction factor
! mu = mu_photoelectric*Fcompton
! where Fcompton=mu_total/mu_photoelectric
! and mu_total is usually mu_photoelectric+mu_compton
! (elastic may be also added)
do j=1,xtal%nbatom
read (25,*,err=79) xtal%FP(j,i),xtal%FPP(j,i),xtal%Fcompton(j,i)
end do
else
do j=1,xtal%nbatom
read (25,*,err=79) xtal%FP(j,i),xtal%FPP(j,i)
xtal%Fcompton(j,i) = 1.0d0
end do
endif
END DO
xtal%filename = FILEIN
CLOSE (25)
RETURN
77 CONTINUE
79 CONTINUE
print *,'crystal_loadCrystalData: Error reading file: '//trim(FILEIN)
stop
END SUBROUTINE crystal_loadCrystalData
!C+++
!C SUBROUTINE crystal_printCrystalData
!C
!C PURPOSE Print the stored crystallographic data from
!C preprocessor (xop/xcrystal or bragg)
!C
!C ALGORITHM
!C
!C MODIFIED Created by M. Sanchez del Rio
!C
!C
!C---
SUBROUTINE crystal_printCrystalData (xtal)
implicit none
type(crystalData),intent(in) :: xtal
integer(kind=ski) :: i,j
print *,' '
print *,'============== in crystal_printCrystalData =================='
print *,'NMAXATOMS=',NMAXATOMS
print *,'NMAXENER=',NMAXENER
print *,'RN=',xtal%rn
print *,'d_spacing=',xtal%d_spacing
print *,'ATNUM=',xtal%atnum(1:xtal%NBATOM)
print *,'TEMPER=',xtal%temper(1:xtal%NBATOM)
print *,'NBATOM,=',xtal%nbatom
print *,'G=',xtal%g(1:xtal%NBATOM)
print *,'G_BAR=',xtal%g_bar(1:xtal%NBATOM)
print *,'G_0=',xtal%g_0(1:xtal%NBATOM)
print *,'NF0COEFF=',xtal%nf0coeff
print *,'F0COEFF=',xtal%f0coeff(1:xtal%NBATOM,1:11)
print *,'NPOINT=',xtal%npoint
print *,'FRACT=',xtal%fract(1:xtal%NBATOM)
print *,'filename=',trim(xtal%filename)
print *,"energy atom_index FP FPP "
print *,'============================================================='
print *,' '
DO I = 1, xtal%NPOINT
do j=1,xtal%nbatom
print *,xtal%energy(i),j, xtal%FP(j,i),xtal%FPP(j,i)
end do
END DO
RETURN
END SUBROUTINE crystal_printCrystalData
!C+++
!C SUBROUTINE CRYSTAL
!C
!C PURPOSE Computes the reflectivity of a symmetric Bragg crystal
!C according to the dynamic theory of x-ray diffraction.
!C
!C ALGORITHM Reference B.E.Warren, X-Ray Diffraction,
!C Addison-Wesley 1969. See also M.J.Bedzyk, G.Materlik
!C and M.V.Kovalchuk, Phy. Rev. B30, 2453(1984).
!C For mosaic crystal reflectivity see Zachariasen,
!C Theory of x-ray diffraction in crystals, Dover (1966)
!C formula 4.24, and Bacon and Lowde, Acta Crystall.1
!C pag 303 (1948) formula 17.
!C
!C MODIFIED July 1989, M. Sanchez del Rio for asymmetry part,
!C July 1990, mosaic part.
!C August 1992, laue part.
!C Dec 2012, cleaned
!C
!C---
!C
!C kwhat: flag:
!C <0 read crystal file (call crystal_loadCrystalData) (note that
!C thio option changes oe1)
!C 0 set mode (call crystal_set)
!C 1 info mode (call crystal_info)
!C 2 reflectivity mode (call crystal_{perfect,mosaic})
SUBROUTINE CRYSTAL (KWHAT, & ! input
oe1, & ! input (also output if kwhat<1)
Q_PHOT, & ! inputs (kwhat>=0)
VIN, VOUT, BH, SURFNOR, & ! inputs (kwhat=2)
R_S, R_P,PHASE_S, PHASE_P, L_EXT_S, L_EXT_P) !, & ! output(kwhat=2)
!THETA_B, SSR, SPR) ! output(kwhat=2)
implicit none
integer(kind=ski),intent(in) :: KWHAT
type(oeSetup),intent(inout) :: oe1
real(kind=skr),intent(in) :: Q_PHOT !Arguments
real(kind=skr),dimension(3),intent(in) :: VIN, VOUT, BH, SURFNOR
real(kind=skr),intent(out) :: R_S, R_P,PHASE_S, PHASE_P
real(kind=skr),intent(out) :: L_EXT_S, L_EXT_P
real(kind=skr) :: ass_fac
real(kind=skr) :: phot,r_lam0
real(kind=skr) :: graze
real(kind=skr) :: absorp
real(kind=skr) :: q_mos
real(kind=skr) :: theta_b_h,theta_b_cor,theta_b_sym
complex(kind=skx):: F_0,FH,FH_BAR,REFRAC,STRUCT
complex(kind=skx):: psi_h,psi_hbar,psi_0
integer(kind=ski):: i_debug=0
! first the usual mode (kwhat=2) to be fast
if (kwhat.eq.2) then
!C***************************************************
!C If flag is = 2, gets reflectivity
!C***************************************************
PHOT = Q_PHOT/TWOPI*TOCM
R_LAM0 = TWOPI/Q_PHOT
!THETA_B = oe1%THETA_B
! these are needed in output
L_EXT_S = oe1%L_EXT_S
L_EXT_P = oe1%L_EXT_P
GRAZE = ASIN(R_LAM0/oe1%crystalData%d_spacing/2.0D0)
if (i_debug.gt.0) write (i_debug,*) '<><>CRYSTAL: working at energy: ',phot
!!C
!!C call crystal_fh to get the structure factor
!!C
if (i_debug.gt.0) write (i_debug,*) '<><>CRYSTAL: calling crystal_fh: ',kwhat
Call crystal_fh (KWHAT, oe1, PHOT, GRAZE , &
FH,FH_BAR,F_0,PSI_H,PSI_HBAR,PSI_0,REFRAC,STRUCT)
if (i_debug.gt.0) then
write(i_debug,*) '<><>CRYSTAL: back crystal_fh, oe1%crystalData%rn: ',oe1%crystalData%rn
write(i_debug,*) 'CRYSTAL: THETA_B=',oe1%THETA_B
write(i_debug,*) 'CRYSTAL: THETA_B[deg]=',oe1%THETA_B*TODEG
write(i_debug,*) 'CRYSTAL: L_EXT_S=',oe1%L_EXT_S
end if
!C
!C Main calculation
!C
IF (oe1%F_MOSAIC.EQ.1) THEN ! mosaic crystal
Q_MOS = pi**2*abs(psi_h*psi_hbar)/r_lam0/sin(2*graze)
!note that crystal)mosaic overwrite the extinction with
!secondary extinction
CALL CRYSTAL_MOSAIC ( Q_PHOT, VIN, BH, SURFNOR, oe1, &
Q_MOS,REFRAC, R_S, R_P,PHASE_S, PHASE_P, L_EXT_S, L_EXT_P)
ELSEif(oe1%F_MOSAIC.EQ.0)then ! perfect crystal
call CRYSTAL_PERFECT (Q_PHOT, VIN, VOUT, BH, SURFNOR,oe1, &
PSI_0,PSI_H,PSI_HBAR, R_S, R_P,PHASE_S, PHASE_P)
ELSEif(oe1%F_MOSAIC.EQ.2)then ! bent, multilamellar
call CRYSTAL_BENT_ML (Q_PHOT, VIN, VOUT, BH, SURFNOR, oe1, &
PSI_0,PSI_H,PSI_HBAR,REFRAC, STRUCT, FH,FH_BAR,R_S, R_P,PHASE_S, PHASE_P)
Elseif(oe1%F_MOSAIC.EQ.3)then ! bent, penning-polder
call CRYSTAL_BENT_PP (Q_PHOT, VIN, VOUT,BH, SURFNOR, oe1, &
PSI_0,PSI_H,PSI_HBAR,STRUCT, REFRAC,R_S, R_P,PHASE_S, PHASE_P)
Elseif(oe1%F_MOSAIC.EQ.4)then
!call bent tt
END IF
return
end if
! then the unusual mode (kwhat < 2)
select case(kwhat)
case(0)
CALL CRYSTAL_SET (Q_PHOT,oe1)
case(1)
CALL CRYSTAL_INFO (oe1)
case default ! load file
Call crystal_loadCrystalData (oe1%FILE_REFL,oe1%crystalData)
if (i_debug.gt.0) Call crystal_printCrystalData (oe1%crystalData)
end select
return
END SUBROUTINE CRYSTAL
!C+++
!C SUBROUTINE CRYSTAL_INFO
!C
!C PURPOSE Prints (on the screen and on fortran unit 23) some
!C info from the crystal set at a given photon energy
!C
!C ALGORITHM
!C
!C MODIFIED M. Sanchez del Rio
!C
!C---
SUBROUTINE CRYSTAL_INFO (oe1)
implicit none
type(oeSetup),intent(inout) :: oe1
real(kind=skr) :: ssr
real(kind=skr) :: spr
real(kind=skr) :: THETA_B
integer(kind=ski):: i
real(kind=skr) :: q_phot,phot,r_lam0
real(kind=skr) :: graze,sin_gra,ass_fac
real(kind=skr) :: absorp
real(kind=skr) :: q_mos
real(kind=skr) :: tmax_mos
real(kind=skr) :: texts,textp
complex(kind=skx):: F_0,FH,FH_BAR,REFRAC,STRUCT
complex(kind=skx):: psi_h,psi_hbar,psi_0
real(kind=skr) :: theta_b_sym,theta_b_h,theta_b_cor
real(kind=skr) :: sigma_gamma,biga_mosaic,omega,kpolsqrt
integer(kind=ski):: i_debug=0
q_phot = oe1%q_phot_set
phot = oe1%phot_set
!C
!C Computes reflectivities at given wavelength and angle.
!C
!! PHOT = Q_PHOT/TWOPI*TOCM
!! R_LAM0 = TWOPI/Q_PHOT
!! SIN_GRA = R_LAM0/oe1%crystalData%d_spacing/2.0D0
!! GRAZE = ASIN(SIN_GRA)
GRAZE = oe1%graze
!todo: be sure that crystal_set has been called before calling crystal_info
!CALL CRYSTAL_SET (Q_PHOT,oe1) !, &
! !THETA_B,THETA_B_H,THETA_B_SYM, THETA_B_COR, L_EXT_S,L_EXT_P,ASS_FAC,SSR,SPR,Q_MOS)
THETA_B = oe1%THETA_B
THETA_B_H = oe1%THETA_B_H
THETA_B_SYM = oe1%THETA_B_SYM
THETA_B_COR = oe1%THETA_B_COR
!L_EXT_S = oe1%L_EXT_S
!L_EXT_P = oe1%L_EXT_P
ASS_FAC = oe1%ASS_FAC
SSR = oe1%SSR
SPR = oe1%SPR
Q_MOS = oe1%Q_MOS
if (i_debug.gt.0) then
write (i_debug,*) '<><>CRYSTAL_INFO ssr: ',ssr
write (i_debug,*) '<><>CRYSTAL_INFO spr: ',spr
write (i_debug,*) '<><>CRYSTAL_INFO theta_b: ',theta_b
end if
!C***************************************************
!C display results on screen (unit=6) and in unit 27
!C***************************************************
!C
!C
!C Calculate mosaic parameters
!C
if (oe1%f_mosaic.EQ.1) then
if (oe1%f_refrac.eq.0) then
biga_mosaic = oe1%thickness*absorp/sin(torad*graze) !bragg
else if (oe1%f_refrac.eq.1) then
biga_mosaic = oe1%thickness*absorp/cos(torad*graze) !laue(alpha=90)
else
write (*,*) 'CRYSTAL_INFO Error: Option not implemented.'
stop
endif
omega = (1/sqrt(twopi))*(1/oe1%spread_mos)
!C
!C Q_mos is the mosaic q variable (Zachariasen's [4.33b])
!C Sigma_Gamma refers to Zachariasen's formula 4.33a
!C
!q_mos = pi**2*abs(psi_h*psi_hbar)/r_lam0/sin(2*graze)
sigma_gamma = omega*q_mos
kpolsqrt = (cos(2*graze))**2
tmax_mos = cos(graze)*dexp(1+2*sigma_gamma/absorp)/2/sigma_gamma/absorp/absorp
!C texts= r_lam0**2/pi/sin(2*graze)/ssr/2/d_spacing
!C textp= r_lam0**2/pi/sin(2*graze)/spr/2/d_spacing
texts= oe1%L_EXT_S*sin(GRAZE)
textp= oe1%L_EXT_P*sin(GRAZE)
end if
DO i=6,23,17 ! write output in units 6(screen) and 23
WRITE(i,*) ' '
call crystal_fh (i, oe1, &
PHOT, GRAZE,FH,FH_BAR,F_0,PSI_H,PSI_HBAR,PSI_0,REFRAC,STRUCT)
WRITE(i,*) 'Photon energy = ',Q_PHOT/TWOPI*TOCM ,' eV'
WRITE(i,*) ' wavelength = ',TWOPI/Q_PHOT*1D8,' Angstroms'
WRITE(i,*) 'Theta (grazing uncorrected Bragg angle) = ',GRAZE,' rads'
WRITE(i,*) ' = ',GRAZE*TODEG,' degrees'
WRITE(i,*) 'Symmetric Bragg angle (corrected) = ', THETA_B_SYM,' rad'
WRITE(i,*) ' = ', THETA_B_SYM*TODEG,' degrees'
WRITE(i,*) 'Angle (H,kin) = Bragg angle (corrected) including asymmetry = ', &
THETA_B_COR, 'rad'
WRITE(i,*) ' = ', &
THETA_B_COR*TODEG,' degrees'
WRITE(i,*) 'Asymmetry angle alpha= ',oe1%a_bragg*todeg,' degrees'
WRITE(i,*) 'Asymmetry factor b= ',ass_fac
WRITE(i,*) ' '
!C Following line modified by T. Jach, 10.10.2001
WRITE(i,*) 'Extinction lengths and depths given here are for amplitude (for intensities, '
WRITE(i,*) ' consider half value).'
WRITE(i,*) 'DEPTH values are measured perpendicular to the crystal surface and LENGTH'
WRITE(i,*) ' values are along the incident beam path.'
WRITE(i,*) 'Extinction length (amplitude) sigma= ', oe1%L_EXT_S*1D4/abs(oe1%vin_bragg(3)),' microns'
!C Following line modified by T. Jach, 10.10.2001
WRITE(i,*) 'Extinction length (amplitude) pi= ', oe1%L_EXT_P*1D4/abs(oe1%vin_bragg(3)),' microns'
!C Following line modified by T. Jach, 10.10.2001
WRITE(i,*) 'Extinction depth (amplitude) sigma= ',oe1%L_EXT_S*1D4,' microns'
!C Following line modified by T. Jach, 10.10.2001
WRITE(i,*) 'Extinction depth (amplitude) pi= ', oe1%L_EXT_P*1D4,' microns'
!C Following line modified by srio 26.10.2001
WRITE(i,*) 'Pendellosung period (amplitude) sigma = Extinction depth sigma times 2 pi = ', oe1%L_EXT_S*1D4*2D0*PI, ' microns'
!C Following line added srio 26.10.2001
WRITE(i,*) 'Pendellosung period (amplitude) pi = Extinction depth pi times 2 pi = ', oe1%L_EXT_P*1D4*2.D0*PI, ' microns'
WRITE(i,*) ' '
if (oe1%f_bragg_a.eq.0) then
WRITE(i,*) 'width of diffraction profile s-pol = ', 2.0D0*SSR*1.0D+6,' microradians'
WRITE(i,*) ' = ', 2.0D0*SSR*TODEG*3600D0,' arc sec'
WRITE(i,*) 'width of diffraction profile p-pol = ', 2.0D0*SPR*1.0D+6,' microradians'
WRITE(i,*) ' = ', 2.0D0*SPR*TODEG*3600D0,' arc sec'
else
WRITE(i,*) 'The width of the diffraction profile as a function of incident angle is'
WRITE(i,*) 'width for s-pol = ', 2.0D0*SSR*1.0D+6,' microradians'
WRITE(i,*) ' = ', 2.0D0*SSR*TODEG*3600D0,' arc sec'
WRITE(i,*) 'width for p-pol = ', 2.0D0*SPR*1.0D+6,' microradians'
WRITE(i,*) ' = ', 2.0D0*SPR*TODEG*3600D0,' arc sec'
!!WRITE(i,*) 'Incident Grazing Angle = ', (oe1%a_bragg+graze)*todeg,' deg'
WRITE(i,*) 'Incident Corrected Angle = ', (theta_b)*todeg,' deg'
WRITE(i,*) ' '
WRITE(i,*) 'The width of the diffraction profile as a function of reflected angle is'
WRITE(i,*) 'width for s-pol = ', 2.0D0*SSR*1.0D+6*(abs(ass_fac)),' microradians'
WRITE(i,*) ' = ', 2.0D0*SSR*(abs(ass_fac))*TODEG*3600D0,' arc sec'
WRITE(i,*) 'width for p-pol = ', 2.0D0*SPR*1.0D+6*(abs(ass_fac)),' microradians'
WRITE(i,*) ' = ', 2.0D0*SPR*(abs(ass_fac))*TODEG*3600D0,' arc sec'
!!WRITE(i,*) 'Reflected Grazing Angle = ', (graze-oe1%a_bragg)*todeg,' deg'
WRITE(i,*) 'Reflected Corrected Angle = ', (theta_b_h)*todeg,' deg'
endif
END DO
i = 23 ! only to file
if (oe1%f_mosaic.EQ.1) then
write (i,*) ' '
write (i,*) '*********** MOSAIC PARAMETERS ***************'
write (i,*) ' '
write (i,*) 'spread_mos= ',2.35*oe1%spread_mos/TORAD ,' deg fwhm'
write (i,*) 'true absorp depth = ',sin(graze)/absorp*1d4, ' microns'
write (i,*) 'true absorp length = ',1d4/absorp,' microns'
write (i,*) 'peak thickness = ',tmax_mos,' cm'
write (i,*) ' '
write (i,*) 'For s-polarization we have: '
write (i,*) ' Q = ',q_mos,' cm^-1 '
!C write (i,*) ' Peak reflectivity = ',r_s
write (i,*) ' Primary extinction:'
write (i,*) ' mu =',sin(graze)/texts,'cm^-1'
write (i,*) ' length =',texts/sin(graze)*1d4,' microns'
write (i,*) ' depth =',texts*1d4,' microns'
write (i,*) ' Secondary extinction:'
write (i,*) ' mu =',sigma_gamma,'cm^-1'
write (i,*) ' length =',1d4/sigma_gamma,' microns'
write (i,*) ' depth =',sin(graze)/sigma_gamma*1d4, ' microns'
write (i,*) ' '
write (i,*) 'For p-polarization we have: '
write (i,*) ' Q = ',q_mos*kpolsqrt,' cm^-1 '
!C write (i,*) ' Peak reflectivity = ',r_p
write (i,*) ' Primary extinction:'
write (i,*) ' mu =',sin(graze)/textp,'cm^-1'
write (i,*) ' length =',textp/sin(graze)*1d4,' microns'
write (i,*) ' depth =',textp*1d4,' microns'
write (i,*) ' Secondary extinction:'
write (i,*) ' mu =',sigma_gamma*kpolsqrt,'cm^-1'
write (i,*) ' length =',1d4/sigma_gamma/kpolsqrt, ' microns'
write (i,*) ' depth =', sin(graze)/sigma_gamma/kpolsqrt*1d4, ' microns'
write (i,*) ' '
endif
write (i,*) ' '
write (i,*) ' '
write (i,*) '*********** VECTOR AND ANGLES ***************'
write (i,*) ' '