-
Notifications
You must be signed in to change notification settings - Fork 88
/
F1F2IN21.f
2698 lines (2050 loc) · 79.6 KB
/
F1F2IN21.f
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
CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
CCC CCC
CCC F1F220 version 0.995 - February 20, 2021 CCC
CCC Collection of subroutines to calculate inclusive cross sections CCC
CCC for range of nuclei. For A > 2 the parameterization is based on CCC
CCC by M. E. Christy, T. Gautam, and A Bodek to 12C, 27Al, 56Fe and CCC
CCC 64Cu. However, the fit scales relatively well with 'A' and CCC
CCC should be good for all nuclei with 10 < A < 80. CCC
CCC Also included is the proton cross section fit and a preliminary CCC
CCC deuteron/neutron fit by M. E. Christy, N. Kalantarians, J. Either CCC
CCC and W. Melnitchouk (to be published) based on both inclusive CCC
CCC deuteron and tagged deuteron data on n/d from BONuS. CCC
CCC Range of validity is W^2 < 32, Q^2 < 32.0 CCC
CCC New data included in deuteron fit, including photoproduction at CCC
CCC Q^2 = 0. CCC
CCC CCC
CCC CCC
CCC A > 2 nuclei fit are 4He, 12C, 27Al, and 56Fe. More nuclei CCC
CCC will be fit in the next version. CCC
CCC CCC
CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
SUBROUTINE F1F2IN21(Z, A, QSQ, WSQ, F1, F2)
!--------------------------------------------------------------------
! Fit to inelastic cross sections for A(e,e')X
! valid for all W^<20 GeV2 and all Q2<30 GeV2
!
! Inputs: Z, A (real*8) are Z and A of nucleus
! (use Z=0., A=1. to get free neutron)
c Qsq (real*8) is 4-vector momentum transfer squared (positive in
c chosen metric)
c Wsq (real*8) is invarinat mass squared of final state calculated
c assuming electron scattered from a free proton
c
c Outputs: F1, F2 (real*8) are structure functions per nucleus
! Version of 02/20/2021 E. Christy
! Made to be consistent with calling F1F2IN09 from P. Bosted
! with much code borrowed.
!--------------------------------------------------------------------
implicit none
real*8 Z,A,QSQ,WSQ
real*8 avgn,F1,F2,FL,F1p,F1n,F2p,F2n,FLp,FLn
real*8 sigm,sigl,sigt,eps
integer IA,IZ,wfn,opt
logical off,doqe/.false./
off = .true.
IA = int(A)
IZ = int(Z)
avgN = A-z
if(IA.LT.2) then
call sf(WSQ,QSQ,F1p,FLp,F2p,F1n,FLn,F2n)
if(IZ.LT.1) then !!! Neutron
F1 = F1n
F2 = F2n
else !!! Proton
F1 = F1p
F2 = F2p
endif
elseif(IA.EQ.2) then !!! Deuteron
eps = 0.5 !!! pass 0 < eps < 1
wfn = 2 !!! CD-Bonn, default
doqe = .false. !!! don't include QE !!!
call rescsd(WSQ,QSQ,eps,doqe,F1,F2,FL,wfn,sigm)
call smearsub(WSQ,QSQ,wfn,off,F2,F1,FL)
c write(6,*) "in f1f2in21: ", wsq,qsq,eps,wfn,f1,f2
elseif(IA.GT.2) then !!! A>2 nuclei
opt = 3 !!! inelastic + MEC
call SFCROSS(WSQ,QSQ,A,Z,opt,sigt,sigl,F1,F2,FL)
endif
return
end
CCC-----------------
SUBROUTINE F1F2QE21(Z, A, QSQ, WSQ, F1qe, F2qe)
!--------------------------------------------------------------------
! Fit to QE cross sections for A(e,e')X
! valid for all W^<30 GeV2 and all Q2<30 GeV2
!
! Inputs: Z, A (real*8) are Z and A of nucleus
! (use Z=0., A=1. to get free neutron)
c Qsq (real*8) is 4-vector momentum transfer squared (positive in
c chosen metric)
c Wsq (real*8) is invarinat mass squared of final state calculated
c assuming electron scattered from a free proton
c
c Outputs: F1, F2 (real*8) are structure functions per nucleus
! Version of 02/20/2021 E. Christy
!
! Note: This is a calling routine in order to be consistent with
! calling of F1F2QE09 from P. Bosted. For the deuteron the
! QE distribution is calculated by smearing with a realistic
! wavefunction in the weak-binding approximation. For A>2
! a superscaling model is used for the fit.
!--------------------------------------------------------------------
IMPLICIT none
real*8 wsq,qsq,A,Z,sigt,sigl,F1qe,F2qe,FLqe
integer opt/1/
integer wfn/2/
logical first/.false./
if(A.EQ.2.0) then
call SQESUB(wsq,qsq,wfn,f2qe,f1qe,fLqe,first)
elseif(A.GT.2.0) then
call sfcross(wsq,qsq,A,Z,opt,sigt,sigl,F1qe,F2qe,FLqe)
endif
c write(6,*) wsq,qsq,a,z,f1qe,f2qe,fLqe
end
CCC-----------------
C=======================================================================
Subroutine FORMFACTS(q2,gmp,gep,gmn,gen)
CCC Returns proton and neutron form factors based on new proton fit CCC
CCC by M.E. Christy including GMp12 data. Neutron starts with CCC
CCC Kelly fit, but includes correction factors extracted from fit to CCC
CCC inclusive deuteron data. Version from July 22, 2020 CCC
!--------------------------------------------------------------------
IMPLICIT NONE
REAL*8 q2,tau,gd,gmp,gep,gmn,gen,mcor,ecor
REAL*8 mu_p/ 2.792782/ ,mu_n/ -1.913148 /
REAL*8 mp/ 0.9382727 /, mp2
mp2 = mp*mp
tau = q2 / 4.0 / mp2
CCC 2021 Christy fit to GMp and GEp including GMp12 data CCC
GMP = mu_p*(1.+0.099481*tau)/
& (1.0+11.089*tau+19.374*tau*tau+5.7798*tau**3)
GEP = (1.+0.24482*tau**2)/
& (1.0+11.715*tau+11.964*tau*tau+27.407*tau**3)
GD = (1./(1 + q2/0.71))**2
CCC 2021 fit to deuteron inclusive data CCC
GMn = mu_n*(1.0+0.12417E-04*tau)/
& (1.000+11.404*tau+17.014*tau*tau+31.219*tau**3)
GEn = (1.5972*tau / (1.0+0.19655*tau)) * GD
return
end
SUBROUTINE RESCSP(w2,q2,sigtp,siglp)
CCCC Returns proton transverse and longitudinal cross sections CCCC
CCCC February 4, 2021 CCCC
IMPLICIT none
real*8 w2,q2,sigtp,siglp
real*8 xvalp(100),xval1(50),xvalL(50)
Integer i
real*8 sigtdis,sigLdis,w2max,w2min
data xvalp /
& 0.12287E+01,0.15196E+01,0.15044E+01,0.17009E+01,0.16765E+01,
& 0.14455E+01,0.12472E+00,0.23000E+00,0.91162E-01,0.88525E-01,
& 0.80876E-01,0.37521E+00,0.76849E+01,0.45411E+01,0.42458E+01,
& 0.18299E+01,0.68593E+01,0.13435E-01,0.17015E+04,0.45076E+01,
& 0.59724E+00,0.18944E+02,0.54924E-01,0.24830E+01,0.30381E+00,
& 0.19419E+02,0.13587E+00,0.20332E+01,0.24378E+01,0.48652E+01,
& 0.30850E+01,0.26876E+01,0.10685E-02,0.97878E+04,0.14193E+00,
& 0.21888E+01,0.25592E+00,0.73812E+00,0.20187E+01,0.93656E-01,
& 0.12867E-01,0.39477E+00,0.19371E+01,0.46900E+01,0.12106E-01,
& 0.17154E+00,0.19929E+01,0.55000E+00,0.42980E+01,0.41895E+00,
& 0.99122E+00,0.99085E+00,0.99798E+00,0.10028E+01,0.98145E+00,
& 0.10163E+01,0.10226E+01,0.10165E+01,0.00000E+00,0.00000E+00,
& 0.00000E+00,0.00000E+00,0.11353E-05,0.10605E+02,0.19705E+01,
& 0.00000E+00,0.19429E-09,0.75543E+02,0.60997E+01,0.00000E+00,
& 0.73839E+01,0.59820E-05,0.19247E+01,0.00000E+00,0.23239E+01,
& 0.16384E+01,0.14220E+01,0.00000E+00,0.10015E-03,0.57993E+00,
& 0.64963E+00,0.00000E+00,0.10093E-03,0.41844E+00,0.42011E+00,
& 0.00000E+00,0.39842E+00,0.86421E+01,0.10000E+00,0.10000E+01,
& -.82435E+01,0.48718E+02,0.51146E+02,0.29711E+01,0.11295E-02,
& 0.47865E+00,0.00000E+00,0.00000E+00,0.00000E+00,0.18935E+02 /
do i=1,50
xval1(i) = xvalp(i)
xvalL(i) = xvalp(50+i)
if(i.LE.12) xvalL(i) = xval1(i)
if(i.EQ.47.OR.i.EQ.48) xvalL(i) = xval1(i)
enddo
w2max = 34.0 !!! max for resonance fit
w2min = 36.0 !!! min for dis fit
call resmodp(1,w2,q2,xval1,sigTp)
call resmodp(2,w2,q2,xvalL,sigLp)
c if(w2.GT.w2min) then
c call disp(w2,q2,sigTdis,sigLdis)
c if(w2.LE.w2max) then
c sigTp = (w2max-w2)*sigTp+(w2-w2min)*sigTdis
c sigLp = (w2max-w2)*sigLp+(w2-w2min)*sigLdis
c sigTp = sigTp/2.0
c sigLp = sigLp/2.0
c else
c sigTp = sigTdis
c sigLp = sigLdis
c endif
c endif
return
end
SUBROUTINE RESCSD(w2,q2,eps,doqe,f1dqe,f2dqe,fLdqe,wfn,sigm)
CCCC Calcultes deuteron cross sections from smeared p+n. CCCC
CCCC Requires SQESUB be called first to read in smearing functions. CCCC
CCCC This should be one at the lowest level to keep from reading CCCC
CCCC multiple times.
IMPLICIT none
real*8 w2,q2,q2t,eps,t1,x,gamma2,q2min
real*8 sigtp,sigtn,siglp,sigln,sigt,sigl,sigm,m
real*8 m2,pi2,alpha,f1d,f2d,fLd,f1dqe,f2dqe,fLdqe
real*8 xvaln(100),off_mKP_fit,delta,xfr,xfl
integer i,j,k,ntssot,wfn,drn
logical doqe,off
c INCLUDE 'parm.cmn'
external off_mKP_fit
off = .false. !! off-shell corrections
c off = .true.
q2min = 4.0E-5
q2t = q2
if(q2.LT.q2min) q2t = q2min
drn = 5
xfr = 0.95
xfl = 1.0E-3
alpha = 1./137.036
m = (0.938272+0.939565)/2.0d0 !!! average p,n
m2 = m*m
pi2 = 3.14158*3.14159
c if(q2.LT.q2min) q2 = q2min !!! Hack to fix photoproduction
x = q2/(w2-m2+q2)
gamma2 = 1.+4.*m2*x*x/q2
t1 = gamma2*f2dqe-2.*x*f1dqe
if(f2dqe.LT.0.0.OR.f1dqe.LT.0.0.OR.fLdqe.LT.0.0)
& write(34,*) w2,q2,f2dqe,f1dqe,fLdqe
call SMEARSUB(w2,q2t,wfn,off,f2d,f1d,fLd)
if(doqe) then
f1d = f1d+f1dqe
f2d = f2d+f2dqe
fLd = fLd+fLdqe
endif
sigt = 0.3894e3*f1d*pi2*alpha*8.0/abs(w2-m2)
sigl = 0.3894e3*fLd*pi2*alpha*8.0/abs(w2-m2)/2.*abs(w2-m2+q2)/q2
sigm = sigt+eps*sigl
return
end
SUBROUTINE RESCSN(w2,q2,sigtn,sigln)
IMPLICIT none
real*8 w2,q2,sigtn,sigln
real*8 xval1(50),xvalL(50)
integer i
real*8 xvaln(100)/
& 0.12287E+01,0.15195E+01,0.15042E+01,0.17062E+01,0.16787E+01,
& 0.14474E+01,0.12502E+00,0.23000E+00,0.90468E-01,0.86480E-01,
& 0.75000E-01,0.38012E+00,0.71001E+01,0.32634E+01,0.10730E+01,
& 0.22650E+01,0.23212E+01,0.18872E+01,0.37882E+01,0.12554E+01,
& 0.29554E+01,0.95547E+01,0.12848E+02,0.27110E+01,0.97586E+03,
& 0.46182E+03,0.25324E+02,0.47504E+05,0.27562E+00,0.20067E+02,
& 0.23779E+00,0.23150E+01,0.28648E+01,0.12912E+04,0.63948E+01,
& 0.61034E+02,0.33488E+00,0.16043E+01,0.24234E+02,0.12297E+00,
& 0.10214E+00,0.52532E+00,0.22671E+01,0.48200E+01,-.77397E-01,
& 0.15751E-01,0.19872E+01,0.54078E+00,0.24385E+01,0.18916E-01,
& 0.10198E+01,0.97431E+00,0.99305E+00,0.98196E+00,0.10425E+01,
& 0.10123E+01,0.98465E+00,0.98416E+00,0.10284E+01,0.98698E+00,
& 0.10162E+01,0.10125E+01,0.11713E+02,0.94676E+02,0.60923E+01,
& 0.46433E-03,0.29789E-01,0.21407E+02,0.29234E+01,0.11574E+02,
& 0.52793E+00,0.41040E+01,0.16980E+01,0.16856E+02,0.35216E+01,
& 0.10865E+02,0.29016E+01,0.31103E+02,0.11927E+02,0.26535E+02,
& 0.81055E+01,0.16015E+01,0.28755E-01,0.24493E+01,0.83532E+00,
& 0.19400E+00,0.15667E+00,0.33667E+01,0.21528E+01,0.10000E+01,
& 0.29325E+02,0.56383E+01,0.73032E+01,0.72695E+01,0.18496E+00,
& 0.78960E+00,0.19850E+01,0.44998E+00,0.10037E+01,0.28875E+01 /
do i=1,50
xval1(i) = xvaln(i)
xvalL(i) = xvaln(50+i)
if(i.LE.12) xvalL(i) = xval1(i) !!! use same masses for L as T !!!
if(i.EQ.47.OR.i.EQ.48) xvalL(i) = xval1(i)
enddo
call resmodn(1,w2,q2,xval1,sigtn)
call resmodn(2,w2,q2,xvalL,sigLn)
return
end
SUBROUTINE RESMODP(sf,w2,q2,xval,sig)
CCC Returns proton transverse (sf=1) and longitudinal (sf=2 Cross sections CCC
CCC Version from February 21, 2021 - Author: M.E. Christy CCC
CCC This routine returns proton photo-absorbtion cross sections CCC
CCC for either transverse or longitudinal photons in units of ub/Sr/Gev. CCC
CCC CCC
CCC Fit form is empirical. Interpret physics from it at your own risk. CCC
CCC replaced 2-pi threshold with eta CCC
IMPLICIT NONE
REAL*8 W,w2,q2,mp,mp2,mpi2,xb,sig,xval(50),mass(7),width(7)
REAL*8 height(7),rescoef(6,4)
REAL*8 nr_coef(3,4),sigr(7),wdif(2),sig_nr,pi2,alpha
REAL*8 mpi,meta,intwidth(7),k,kcm,kr(7),kcmr(7),ppicm,ppi2cm
REAL*8 petacm,ppicmr(7),ppi2cmr(7),petacmr(7),epicmr(7),epi2cmr(7)
REAL*8 eetacmr(7),epicm,epi2cm,eetacm,br(7,3),ang(7)
REAL*8 pgam(7),pwid(7,3),x0(7),dip,mon,q20,A0
REAL*8 sig_res,xpr(2),t1,t2
INTEGER i,j,num,sf
mp = 0.9382727
mpi = 0.134977
mpi2 = mpi*mpi
meta = 0.547862
mp2 = mp*mp
alpha = 1./137.036
W = sqrt(w2)
wdif(1) = w - (mp + mpi)
wdif(2) = w - (mp + meta)
q20 = 0.05
q20= xval(50)
CCCC single pion branching ratios CCCC
br(1,1) = 1.00 !!! P33(1232)
br(2,1) = 0.45 !!! S11(1535)
br(3,1) = 0.60 !!! D13(1520)
br(4,1) = 0.65 !!! F15(1680)
br(5,1) = 0.60 !!! S11(1650)
br(6,1) = 0.65 !!! P11(1440) roper
br(7,1) = 0.60 !!! F37(1950)
CCCC eta branching ratios CCCC
br(1,3) = 0.0 !!! P33(1232)
br(2,3) = 0.40 !!! S11(1535)
br(3,3) = 0.08 !!! D13(1520)
br(4,3) = 0.0 !!! F15(1680)
br(5,3) = 0.20 !!! S11(1650)
br(6,3) = 0.0 !!! P11(1440) roper
br(7,3) = 0.0 !!! F37(1950)
CCCC 2-pion branching ratios CCCC
do i=1,7
br(i,2) = 1.-br(i,1)-br(i,3)
enddo
CCCC Meson angular momentum CCCC
ang(1) = 1. !!! P33(1232)
ang(2) = 0. !!! S11(1535)
ang(3) = 2. !!! D13(1520)
ang(4) = 3. !!! F15(1680)
ang(5) = 0. !!! S15(1650)
ang(6) = 1. !!! P11(1440) roper
ang(7) = 3. !!! F37(1950)
do i=1,7 !!! resonance damping parameter !!!
x0(i) = 0.178 !!!
enddo
c x0(1) = 0.125
x0(1) = 0.142
do i=1,7
br(i,2) = 1.-br(i,1)-br(i,3)
enddo
dip = 1./(1.+q2/1.05)**2. !!! Dipole parameterization !!!
mon = 1./(1.+q2/1.05)**1.
xb = q2/(q2+w2-mp2)
xpr(1) = 1.+(w2-(mp+mpi)**2)/(q2+q20)
xpr(1) = 1./xpr(1)
xpr(2) = 1.+(w2-(mp+meta)**2)/(q2+q20)
xpr(2) = 1./xpr(2)
if(w.LE.(mp+mpi)) xpr(1) = 1.0
if(w.LE.(mp+meta)) xpr(2) = 1.0
CCC Calculate kinematics needed for threshold Relativistic B-W CCC
k = (w2 - mp2)/2./mp
kcm = (w2-mp2)/2./w
epicm = (W2 + mpi**2 -mp2 )/2./w
ppicm = SQRT(MAX(0.0,(epicm**2 - mpi**2)))
epi2cm = (W2 + (2.*mpi)**2 -mp2 )/2./w
ppi2cm = SQRT(MAX(0.0,(epi2cm**2 - (2.*mpi)**2)))
eetacm = (W2 + meta*meta -mp2 )/2./w
petacm = SQRT(MAX(0.0,(eetacm**2 - meta**2)))
num = 0
do i=1,6 !!! Read in resonance masses !!!
num = num + 1
mass(i) = xval(i)
enddo
do i=1,6 !!! Read in resonance widths !!!
num = num + 1
intwidth(i) = xval(num)
width(i) = intwidth(i)
enddo
mass(7) = xval(47)
intwidth(7) = xval(48)
width(7) = intwidth(7)
do i=1,7
kr(i) = (mass(i)**2-mp2)/2./mp
kcmr(i) = (mass(i)**2.-mp2)/2./mass(i)
epicmr(i) = (mass(i)**2 + mpi**2 -mp2 )/2./mass(i)
ppicmr(i) = SQRT(MAX(0.0,(epicmr(i)**2 - mpi**2)))
epi2cmr(i) = (mass(i)**2 + (2.*mpi)**2 -mp2 )/2./mass(i)
ppi2cmr(i) = SQRT(MAX(0.0,(epi2cmr(i)**2 - (2.*mpi)**2)))
eetacmr(i) = (mass(i)**2 + meta*meta -mp2 )/2./mass(i)
petacmr(i) = SQRT(MAX(0.0,(eetacmr(i)**2 - meta**2)))
CCC Calculate partial widths CCC
pwid(i,1) = intwidth(i)*(ppicm/ppicmr(i))**(2.*ang(i)+1.)
& *((ppicmr(i)**2+x0(i)**2)/(ppicm**2+x0(i)**2))**ang(i)
c !!! 1-pion decay mode
pwid(i,2) = intwidth(i)*(ppi2cm/ppi2cmr(i))**(2.*ang(i)+4.)
& *((ppi2cmr(i)**2+x0(i)**2)/(ppi2cm**2+x0(i)**2))
& **(ang(i)+2) !!! 2-pion decay mode
pwid(i,2) = W/mass(i)*pwid(i,2)
pwid(i,3) = 0. !!! eta decay mode
if(i.EQ.2.OR.i.EQ.5) then
pwid(i,3) = intwidth(i)*(petacm/petacmr(i))**(2.*ang(i)+1.)
& *((petacmr(i)**2+x0(i)**2)/(petacm**2+x0(i)**2))**ang(i)
c !!! eta decay only for S11's
endif
pgam(i) = (kcm/kcmr(i))**2*
& (kcmr(i)**2+x0(i)**2)/(kcm**2+x0(i)**2)
pgam(i) = intwidth(i)*pgam(i)
width(i) = br(i,1)*pwid(i,1)+br(i,2)*pwid(i,2)+br(i,3)*pwid(i,3)
enddo
CCC End resonance kinematics and Widths calculations CCC
CCC Begin resonance Q^2 dependence calculations CCC
do i=1,6
do j=1,4
num = num + 1
rescoef(i,j)=xval(num)
enddo
if(sf.EQ.1) then
height(i) = rescoef(i,1)*
& (1.+rescoef(i,2)*q2/(1.+rescoef(i,3)*q2))*
& mon**rescoef(i,4)
else
height(i) = (rescoef(i,1)+rescoef(i,2)*q2)
& *exp(-1.*rescoef(i,3)*q2)
endif
height(i) = height(i)*height(i)
enddo
if(sf.EQ.2) then !!! 4th resonance region !!!
height(7) = (xval(44)+xval(45)*q2)*exp(-1.0*xval(46)*q2)
else
height(7) = xval(49)*mon
endif
height(7) = height(7)*height(7)
CCC End resonance Q^2 dependence calculations CCC
do i=1,3 !!! Non-Res coefficients !!!
do j=1,4
num = num + 1
nr_coef(i,j)=xval(num)
enddo
enddo
CCC Calculate Breit-Wigners for all resonances CCC
sig_res = 0.0
do i=1,7
sigr(i) = width(i)*pgam(i)/((W2 - mass(i)**2.)**2.
& + (mass(i)*width(i))**2.)
sigr(i) = height(i)*kr(i)/k*kcmr(i)/kcm*sigr(i)/intwidth(i)
sig_res = sig_res + sigr(i)
enddo
sig_res = sig_res*w
if(sf.EQ.2) sig_res = sig_res*q2
CCC Finish resonances / start non-res background calculation CCC
sig_nr = 0.
if(sf.EQ.1.and.xpr(1).LT.1.0) then
A0 = (1.0+xval(44)*q2)/(1.+q2/xval(42))**xval(43) !!! overall amplitude
t1 = xval(38)*log(2.0+q2/xval(39)) !!! exponent of (1-xpr)
t2 = xval(40)*log(2.0+q2/xval(41))+xval(45) !!! exponent of xpr
if(xpr(1).LE.1.0) then !!! 1-pi threshold
sig_nr = xval(37)*389.4*A0*(1.-xpr(1))**t1*xpr(1)**t2
endif
if(xpr(2).LE.1.0) then !!! eta threshold
sig_nr = sig_nr+xval(46)*389.4*A0*(1.-xpr(2))**t1*xpr(2)**t2
endif
elseif(sf.EQ.2.and.xpr(1).LT.1.0) then
t1 = xval(38)*log(1.001+q2/q20) !!! exponent of (1-xpr)
t2 = xval(41)/(1.001+q2/xval(42))**xval(43) !!! exponent of xpr
if(xpr(1).LE.1.0) then
sig_nr = sig_nr + xval(37)*389.4*
& xb*(1.-xpr(1))**t1*xpr(1)**t2
endif
sig_nr =sig_nr/log(1.001+q2/xval(39))
endif
sig = sig_res + sig_nr
if((w-mp).LT.wdif(1)) sig = 0.0
1000 format(8f12.5)
1001 format(7f12.3)
RETURN
END
SUBROUTINE RESMODN(sf,w2,q2,xval,sig)
CCC Returns proton transverse (sf=1) and longitudinal (sf=2 Cross sections CCC
CCC Version from February 21, 2021 - Author: M.E. Christy CCC
CCC This routine returns proton photo-absorbtion cross sections CCC
CCC for either transverse or longitudinal photons in units of ub/Sr/Gev. CCC
CCC CCC
CCC Fit form is empirical. Interpret physics from it at your own risk. CCC
CCC replaced 2-pi threshold with eta CCC
IMPLICIT NONE
REAL*8 W,w2,q2,mp,mp2,mpi2,xb,sig,xval(50),mass(7),width(7)
REAL*8 height(7),rescoef(6,4)
REAL*8 nr_coef(3,4),sigr(7),wdif(2),sig_nr,pi2,alpha
REAL*8 mpi,meta,intwidth(7),k,kcm,kr(7),kcmr(7),ppicm,ppi2cm
REAL*8 petacm,ppicmr(7),ppi2cmr(7),petacmr(7),epicmr(7),epi2cmr(7)
REAL*8 eetacmr(7),epicm,epi2cm,eetacm,br(7,3),ang(7)
REAL*8 pgam(7),pwid(7,3),x0(7),dip,mon,q20,A0
REAL*8 sig_res,xpr(2),t1,t2
INTEGER i,j,num,sf
mp = 0.939565
mpi = 0.134977
mpi2 = mpi*mpi
meta = 0.547862
mp2 = mp*mp
alpha = 1./137.036
W = sqrt(w2)
wdif(1) = w - (mp + mpi)
wdif(2) = w - (mp + meta)
q20 = 0.05
q20= xval(50)
CCCC single pion branching ratios CCCC
br(1,1) = 1.00 !!! P33(1232)
br(2,1) = 0.45 !!! S11(1535)
br(3,1) = 0.60 !!! D13(1520)
br(4,1) = 0.65 !!! F15(1680)
br(5,1) = 0.60 !!! S11(1650)
br(6,1) = 0.65 !!! P11(1440) roper
br(7,1) = 0.60 !!! F37(1950)
CCCC eta branching ratios CCCC
br(1,3) = 0.0 !!! P33(1232)
br(2,3) = 0.40 !!! S11(1535)
br(3,3) = 0.08 !!! D13(1520)
br(4,3) = 0.0 !!! F15(1680)
br(5,3) = 0.20 !!! S11(1650)
br(6,3) = 0.0 !!! P11(1440) roper
br(7,3) = 0.0 !!! F37(1950)
CCCC 2-pion branching ratios CCCC
do i=1,7
br(i,2) = 1.-br(i,1)-br(i,3)
enddo
CCCC Meson angular momentum CCCC
ang(1) = 1. !!! P33(1232)
ang(2) = 0. !!! S11(1535)
ang(3) = 2. !!! D13(1520)
ang(4) = 3. !!! F15(1680)
ang(5) = 0. !!! S15(1650)
ang(6) = 1. !!! P11(1440) roper
ang(7) = 3. !!! F37(1950)
do i=1,7 !!! resonance damping parameter !!!
x0(i) = 0.178 !!!
enddo
c x0(1) = 0.125
x0(1) = 0.142
do i=1,7
br(i,2) = 1.-br(i,1)-br(i,3)
enddo
dip = 1./(1.+q2/1.05)**2. !!! Dipole parameterization !!!
mon = 1./(1.+q2/1.1)**1.
xb = q2/(q2+w2-mp2)
xpr(1) = 1.+(w2-(mp+mpi)**2)/(q2+q20)
xpr(1) = 1./xpr(1)
xpr(2) = 1.+(w2-(mp+meta)**2)/(q2+q20)
xpr(2) = 1./xpr(2)
if(w.LE.(mp+mpi)) xpr(1) = 1.0
if(w.LE.(mp+meta)) xpr(2) = 1.0
CCC Calculate kinematics needed for threshold Relativistic B-W CCC
k = (w2 - mp2)/2./mp
kcm = (w2-mp2)/2./w
epicm = (W2 + mpi**2 -mp2 )/2./w
ppicm = SQRT(MAX(0.0,(epicm**2 - mpi**2)))
epi2cm = (W2 + (2.*mpi)**2 -mp2 )/2./w
ppi2cm = SQRT(MAX(0.0,(epi2cm**2 - (2.*mpi)**2)))
eetacm = (W2 + meta*meta -mp2 )/2./w
petacm = SQRT(MAX(0.0,(eetacm**2 - meta**2)))
num = 0
do i=1,6 !!! Read in resonance masses !!!
num = num + 1
mass(i) = xval(i)
enddo
do i=1,6 !!! Read in resonance widths !!!
num = num + 1
intwidth(i) = xval(num)
width(i) = intwidth(i)
enddo
mass(7) = xval(47)
intwidth(7) = xval(48)
width(7) = intwidth(7)
do i=1,7
kr(i) = (mass(i)**2-mp2)/2./mp
kcmr(i) = (mass(i)**2.-mp2)/2./mass(i)
epicmr(i) = (mass(i)**2 + mpi**2 -mp2 )/2./mass(i)
ppicmr(i) = SQRT(MAX(0.0,(epicmr(i)**2 - mpi**2)))
epi2cmr(i) = (mass(i)**2 + (2.*mpi)**2 -mp2 )/2./mass(i)
ppi2cmr(i) = SQRT(MAX(0.0,(epi2cmr(i)**2 - (2.*mpi)**2)))
eetacmr(i) = (mass(i)**2 + meta*meta -mp2 )/2./mass(i)
petacmr(i) = SQRT(MAX(0.0,(eetacmr(i)**2 - meta**2)))
CCC Calculate partial widths CCC
pwid(i,1) = intwidth(i)*(ppicm/ppicmr(i))**(2.*ang(i)+1.)
& *((ppicmr(i)**2+x0(i)**2)/(ppicm**2+x0(i)**2))**ang(i)
c !!! 1-pion decay mode
pwid(i,2) = intwidth(i)*(ppi2cm/ppi2cmr(i))**(2.*ang(i)+4.)
& *((ppi2cmr(i)**2+x0(i)**2)/(ppi2cm**2+x0(i)**2))
& **(ang(i)+2) !!! 2-pion decay mode
pwid(i,2) = W/mass(i)*pwid(i,2)
pwid(i,3) = 0. !!! eta decay mode
if(i.EQ.2.OR.i.EQ.5) then
pwid(i,3) = intwidth(i)*(petacm/petacmr(i))**(2.*ang(i)+1.)
& *((petacmr(i)**2+x0(i)**2)/(petacm**2+x0(i)**2))**ang(i)
c !!! eta decay only for S11's
endif
pgam(i) = (kcm/kcmr(i))**2*
& (kcmr(i)**2+x0(i)**2)/(kcm**2+x0(i)**2)
pgam(i) = intwidth(i)*pgam(i)
width(i) = br(i,1)*pwid(i,1)+br(i,2)*pwid(i,2)+br(i,3)*pwid(i,3)
enddo
CCC End resonance kinematics and Widths calculations CCC
CCC Begin resonance Q^2 dependence calculations CCC
do i=1,6
do j=1,4
num = num + 1
rescoef(i,j)=xval(num)
enddo
if(sf.EQ.1) then
height(i) = rescoef(i,1)*
& (1.+rescoef(i,2)*q2/(1.+rescoef(i,3)*q2))*
& mon**rescoef(i,4)
else
height(i) = (rescoef(i,1)+rescoef(i,2)*q2)
& *exp(-1.*rescoef(i,3)*q2)
endif
height(i) = height(i)*height(i)
enddo
if(sf.EQ.2) then !!! 4th resonance region !!!
height(7) = (xval(44)+xval(45)*q2)*exp(-1.0*xval(46)*q2)
else
height(7) = xval(49)*mon
endif
height(7) = height(7)*height(7)
CCC End resonance Q^2 dependence calculations CCC
do i=1,3 !!! Non-Res coefficients !!!
do j=1,4
num = num + 1
nr_coef(i,j)=xval(num)
enddo
enddo
CCC Calculate Breit-Wigners for all resonances CCC
sig_res = 0.0
do i=1,7
sigr(i) = width(i)*pgam(i)/((W2 - mass(i)**2.)**2.
& + (mass(i)*width(i))**2.)
sigr(i) = height(i)*kr(i)/k*kcmr(i)/kcm*sigr(i)/intwidth(i)
sig_res = sig_res + sigr(i)
enddo
sig_res = sig_res*w
if(sf.EQ.2) sig_res = sig_res*q2
CCC Finish resonances / start non-res background calculation CCC
sig_nr = 0.
if(sf.EQ.1.and.xpr(1).LT.1.0) then
A0 = (1.0+xval(44)*q2)/(1.+q2/xval(42))**xval(43) !!! overall amplitude
t1 = xval(38)*log(2.51+q2/xval(39)) !!! exponent of (1-xpr)
t2 = xval(40)*log(2.51+q2/xval(41))+xval(45) !!! exponent of xpr
if(xpr(1).LE.1.0) then !!! 1-pi threshold
sig_nr = xval(37)*389.4*A0*(1.-xpr(1))**t1*xpr(1)**t2
endif
if(xpr(2).LE.1.0) then !!! eta threshold
sig_nr = sig_nr+xval(46)*389.4*A0*(1.-xpr(2))**t1*xpr(2)**t2
endif
elseif(sf.EQ.2.and.xpr(1).LT.1.0) then
t1 = xval(38)*log(1.001+q2/q20) !!! exponent of (1-xpr)
t2 = xval(41)/(1.001+q2/xval(42))**xval(43) !!! exponent of xpr
if(xpr(1).LE.1.0) then
sig_nr = sig_nr + xval(37)*389.4*
& xb*(1.-xpr(1))**t1*xpr(1)**t2
endif
sig_nr =sig_nr/log(1.001+q2/xval(39))
endif
sig = sig_res + sig_nr
if((w-mp).LT.wdif(1)) sig = 0.0
1000 format(8f12.5)
1001 format(7f12.3)
RETURN
END
SUBROUTINE SMEARSUB(w2,q2,wfn,off,f2s,f1s,fLs)
CCCC Fermi smears structure functions. Requires subroutne GETFY CCCC
CCCC to be initialized first. CCCC
IMPLICIT none
real*8 w2,q2,mp,mp2,x,z,wwz,alpha,alpha2,gamma,gamma2
real*8 y,fy11,fy12,fy2,inc,f2s,f1s,fLs,f2pi,f2ni,xvaln(100)
real*8 sigtp,sigtn,siglp,sigln,f1pi,f1ni,fLpi,fLni
real*8 f1i(1000),f2i(1000),vfact1(1000),vfact2(1000)
real*8 vfact1off(1000),vfact2off(1000),fy11off,fy12off
real*8 fy2off,f1soff,f2soff,fLsoff
real*8 ymin,ymax,pi,pi2,epsd
real*8 c0,c1,c2,delta(1000)
integer i,j,nbins,wfn
c logical first/.true./
logical firsty/.false./
logical off
real*8 dsimps
external dsimps
C *****Off-shell parameters*****
if(wfn.eq.1) then
c0=-3.6735d0
c1=0.057717d0
c2=0.36419d0
elseif(wfn.eq.2) then
c0=-7.2061d0
c1=0.062022d0
c2=0.38018d0
elseif(wfn.eq.3) then
c0=1.7204d0
c1=0.050923d0
c2=0.34360d0
elseif(wfn.eq.4) then
c0=-1.7690d0
c1=0.058321d0
c2=0.36976d0
else
c0=0.d0
c1=0.d0
c2=0.d0
write(6,*) "Warning: wavefunction undefined"
endif
C ******************************
mp = (0.938272+0.939565)/2.0d0 !!! average of p, n
mp2 = mp*mp
alpha = 1/137.03599
alpha2 = alpha*alpha
pi = 3.14159
pi2 = pi*pi
epsd = -2.2E-03
x = q2/(w2-mp2+q2)
f2s = 0.0
f1s = 0.0
fLs = 0.0
if(x.LT.0.0D0) return
gamma2 = (1.+4.*mp2*x*x/q2)
gamma = sqrt(gamma2)
nbins = 60
if(w2.LT.1.6) nbins = 80
c nbins = 200 !!! adjust as needed
if(w2.GE.2.6.AND.Q2.GE.2.0) nbins = 46
CCC/// Fill array for Simpson's rule calculation ///CCC
c if(x.LT.1.0) then
c ymin = max(x*(1.0D0-2.0D0*epsd*mp/(w2-mp2)),x)
c else
c ymin = x
c endif
CCC Note problem if Q2 = 0 CCC
ymin = max(x*(1.0D0-2.0D0*epsd*mp/q2),x)
c ymin = min(ymin,2.0d0)
ymax = min(1.0d0+gamma2/2.0d0+epsd/mp,2.0d0)
c write(6,*) "W2, ymin, ymax: ",w2,ymin,ymax
if(ymin.GE.ymax) return
c ymax = 2.0d0
c write(6,*) off
if(ymin.EQ.0) write(6,*) "Error ymin = 0"
inc = (ymax-ymin)/float(nbins)
y = ymin-inc
!$OMP PARALLEL DO PRIVATE(i,y,z,wwz,F1pi,FLpi,f1ni,FLni,f2pi,f2ni,fy11,fy12,fy2)
do i=1,nbins+1 !!!! First calculate the cross smearing function for each y bin !!!!
y = y + inc
z = x/y