forked from JeffersonLab/simc_gfortran
-
Notifications
You must be signed in to change notification settings - Fork 5
/
F1F209.f
1684 lines (1453 loc) · 62.3 KB
/
F1F209.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
ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
c F1F209.f
c Package of FORTRAN subroutine describing fits to inclusive inelastic
c electron scattering from proton, neutron, deuteron, or heaviere nuclei
c Proton fit is described in:
c M.E. Christy and P.E. Bosted, ``Empirical Fit to Precision
c Inclusive Electron-Proton Cross Sections in the Resonance Region'',
c (arXiv:0712.3731). Submitted to Phys. Rev. C.
c Deuteron (and netron) fit is described in:
c P.E. Bosted and M.E. Christy, ``Empirical Fit to Inelastic
c Electron-Deuteron and Electron-Neutron
c Resonance Region Transverse Cross Sections,
c (arXiv:0711.0159), publichished in Phys. Rev. C 77, 065206 (2008). (
c New fits for A>2 by Vahe M. and Peter B. (to be publsihed).
c
c routine in this package:
c Main program test. Compiling and linking all code
c in this package will make a "test" program.
c Upon execution, the output should be:
c 1.000 1.500 0.617 0.204 0.231 0.312 0.410 1.644 0.072 0.102
c 1.000 2.000 0.472 0.129 0.159 0.292 0.250 2.329 0.199 0.191
c 1.000 2.500 0.382 0.169 0.227 0.417 0.360 2.973 0.241 0.292
c 1.000 3.000 0.321 0.224 0.315 0.511 0.500 3.562 0.327 0.359
c 1.000 3.500 0.276 0.218 0.309 0.517 0.521 4.060 0.348 0.387
c 2.000 1.500 0.763 0.079 0.091 0.110 0.104 0.400 0.154 0.188
c 2.000 2.000 0.641 0.058 0.091 0.157 0.104 0.764 0.184 0.216
c 2.000 2.500 0.553 0.090 0.139 0.241 0.162 1.174 0.215 0.246
c 2.000 3.000 0.485 0.137 0.214 0.324 0.252 1.558 0.233 0.270
c 2.000 3.500 0.433 0.136 0.227 0.354 0.272 1.919 0.267 0.283
c 1.000 0.747 0.432
c 2.000 0.157 0.132
c 3.000 0.048 0.049
c
c
c subroutine F1F2IN09. Returns inelastic F1, F2, and R for
c nucleus with charge Z and atomic number A
c for given value of Q2 and W**2. F1 and F2
c are per nucleus (not per nucleon).
c subroutine pind. Returns F1, R, Sigma_L, and Sigma_T
c for a proton with Fermi motion of a detueron
c subroutine resd. Returns F1 for a deuteron.
c subroutine resder. Returns error on F1 from fit.
c Requires auxillary file "F1F207D2emat.dat"
c subroutine resmodd. Returns F1 for average of a free
c proton and neutron
c subroutine christy507. Returns F1, R, sigma_L, and
c sigma_T for the free proton.
c subroutine resmod507. Returns sigma_L or sigma_T
c (called by christy507)
c subroutine mec. Called by resd to get extra terms
c for dip region between quasi-elastic and Delta
c funtion fitemc. Fit to "EMC effect" used to
c get F1 and F2 for A>2
c
c subroutine F1F2QE09. Returns quasi-elastic F1, F2 for
c nucleus with charge Z and atomic number A
c for given value of Q2 and W**2
c
cc implicit none
cc integer iq,iw
cc real*8 q2,w2,F1n,F2n,r,F1p,F2p,F1d,F2d,F1c,F2c
cc real*8 f1dqe, f2dqe, nu, x, am/0.9383/, F1be, F2be,rbe
cc
cc do iq=1,2
cc q2 = 1.0 * iq
cc do iw=1,5
cc w2 = 1. + 0.5*iw
cc nu = (w2 - am**2 + q2) / 2. / am
cc x = q2 / 2. / am / nu
cc call F1F2IN09(0.D0, 1.D0, q2, w2, F1n, F2n,r)
cc call F1F2IN09(1.D0, 1.D0, q2, w2, F1p, F2p,r)
cc call F1F2IN09(1.D0, 2.D0, q2, w2, F1d, F2d,r)
cc call F1F2IN09(4.D0, 9.D0, q2, w2, F1be, F2be,rbe)
cc write(6,'(10f7.3)') q2,w2,x,f2n,f2p,f2d,
cc > f1p,f1be,r,rbe
cc enddo
cc enddo
cc do iq=1,3
cc q2 = 1.0 * iq
cc w2 = am**2
cc call F1F2QE09(1.D0, 2.D0, q2, w2, F1dqe, F2dqe)
cc write(6,'(3f8.3)') q2, F1dqe, F2dqe
cc enddo
cc stop
cc end
C=======================================================================
SUBROUTINE F1F2IN09(Z, A, QSQ, Wsq, F1, F2)
!--------------------------------------------------------------------
! Fit to inelastic cross sections for A(e,e')X
! valid for all W<3 GeV and all Q2<10 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 10/20/2006 P. Bosted
!--------------------------------------------------------------------
implicit none
REAL*8 P(0:23)
COMMON/PARCORR/P
real*8 Z,A,qsq,wsq,w,f1c,x
real*8 avgn,r,dr,nu,eps,kappa,sigres,flux,siglp,sigtp,F1pp,F1dp
REAL*8 W1,W2,sigt,rc,w1pb,w2pb,F1,F2,sigl,F1d, F1p,qv
REAL*8 W1p,W2P,DW2DPF,WSQP,pf,kf,es,dw2des,fyuse
real a4, x4, fitemc, emcfac
logical goodfit
INTEGER ISM
! new variables for Fermi smearing over +/- 3 sigma. Sum of
! fy values is 1.000, so no effect if W1 is constant with W
REAL*8 XX(15)/-3.000,-2.571,-2.143,-1.714,-1.286,-0.857,
> -0.429, 0.000, 0.429, 0.857, 1.286, 1.714,
> 2.143, 2.571, 3.000/
real*8 fy(15)/0.0019, 0.0063, 0.0172, 0.0394, 0.0749, 0.1186,
> 0.1562, 0.1712, 0.1562, 0.1186, 0.0749, 0.0394,
> 0.0172, 0.0063, 0.0019/
! This is for exp(-xx**2/2.), from teste.f
real*8 xxp(99)/
> -3.000,-2.939,-2.878,-2.816,-2.755,-2.694,-2.633,-2.571,-2.510,
> -2.449,-2.388,-2.327,-2.265,-2.204,-2.143,-2.082,-2.020,-1.959,
> -1.898,-1.837,-1.776,-1.714,-1.653,-1.592,-1.531,-1.469,-1.408,
> -1.347,-1.286,-1.224,-1.163,-1.102,-1.041,-0.980,-0.918,-0.857,
> -0.796,-0.735,-0.673,-0.612,-0.551,-0.490,-0.429,-0.367,-0.306,
> -0.245,-0.184,-0.122,-0.061, 0.000, 0.061, 0.122, 0.184, 0.245,
> 0.306, 0.367, 0.429, 0.490, 0.551, 0.612, 0.673, 0.735, 0.796,
> 0.857, 0.918, 0.980, 1.041, 1.102, 1.163, 1.224, 1.286, 1.347,
> 1.408, 1.469, 1.531, 1.592, 1.653, 1.714, 1.776, 1.837, 1.898,
> 1.959, 2.020, 2.082, 2.143, 2.204, 2.265, 2.327, 2.388, 2.449,
> 2.510, 2.571, 2.633, 2.694, 2.755, 2.816, 2.878, 2.939, 3.000/
! these are 100x bigger for convenience
real*8 fyp(99)/
> 0.0272,0.0326,0.0390,0.0464,0.0551,0.0651,0.0766,0.0898,0.1049,
> 0.1221,0.1416,0.1636,0.1883,0.2159,0.2466,0.2807,0.3182,0.3595,
> 0.4045,0.4535,0.5066,0.5637,0.6249,0.6901,0.7593,0.8324,0.9090,
> 0.9890,1.0720,1.1577,1.2454,1.3349,1.4254,1.5163,1.6070,1.6968,
> 1.7849,1.8705,1.9529,2.0313,2.1049,2.1731,2.2350,2.2901,2.3379,
> 2.3776,2.4090,2.4317,2.4454,2.4500,2.4454,2.4317,2.4090,2.3776,
> 2.3379,2.2901,2.2350,2.1731,2.1049,2.0313,1.9529,1.8705,1.7849,
> 1.6968,1.6070,1.5163,1.4254,1.3349,1.2454,1.1577,1.0720,0.9890,
> 0.9090,0.8324,0.7593,0.6901,0.6249,0.5637,0.5066,0.4535,0.4045,
> 0.3595,0.3182,0.2807,0.2466,0.2159,0.1883,0.1636,0.1416,0.1221,
> 0.1049,0.0898,0.0766,0.0651,0.0551,0.0464,0.0390,0.0326,0.0272/
integer iz,ia,i
real PM/0.93828/,THCNST/0.01745329/,ALPHA/137.0388/
real*8 PI/3.1415926535D0/
! deuteron fit parameters
real*8 xvald0(50)/
> 0.1964E+01, 0.1086E+01, 0.5313E-02, 0.1265E+01, 0.8000E+01,
> 0.2979E+00, 0.1354E+00, 0.2200E+00, 0.8296E-01, 0.9578E-01,
> 0.1094E+00, 0.3794E+00, 0.8122E+01, 0.5189E+01, 0.3290E+01,
> 0.1870E+01, 0.6110E+01,-0.3464E+02, 0.9000E+03, 0.1717E+01,
> 0.4335E-01, 0.1915E+03, 0.2232E+00, 0.2119E+01, 0.2088E+01,
> -0.3029E+00, 0.2012E+00, 0.1104E-02, 0.2276E-01,-0.4562E+00,
> 0.2397E+00, 0.1204E+01, 0.2321E-01, 0.5419E+03, 0.2247E+00,
> 0.2168E+01, 0.2266E+03, 0.7649E-01, 0.1457E+01, 0.1318E+00,
> -0.7534E+02, 0.1776E+00, 0.1636E+01, 0.1350E+00,-0.5596E-02,
> 0.5883E-02, 0.1934E+01, 0.3800E+00, 0.3319E+01, 0.1446E+00/
cc
real*8 F1M
logical DEBUG/.TRUE./
IA = int(A)
avgN = A - Z
nu = (wsq - pm**2 + qsq) / 2. / pm
qv = sqrt(nu**2 + qsq)
if(Wsq.le.0.0) W = 0.0
W = sqrt(Wsq)
x = QSQ / (2.0 * pm * nu)
if(Wsq.le.0.0) x = 0.0
! Cross section for proton or neutron
W1 = 0.
W2 = 0.
IF(IA .lt. 2 .and. wsq.gt.1.155) THEN
call CHRISTY507(Wsq,Qsq,F1p,Rc,sigt,sigl)
! If neutron, subtract proton from deuteron. Factor of two to
! convert from per nucleon to per deuteron
if(Z .lt. 0.5) then
call resmodd(wsq,qsq,xvald0,F1d)
F1p = F1d * 2.0 - F1p
endif
W1 = F1p / PM
W2 = W1 * (1. + Rc) / (1. + nu**2 / qsq)
ENDIF
! For deuteron
if(IA .eq. 2) then
c get Fermi-smeared R from Erics proton fit
call pind(Wsq, Qsq, F1c, Rc, sigt, sigl)
c get fit to F1 in deuteron, per nucleon
call resd(qsq, wsq, xvald0, F1d)
c convert to W1 per deuteron
W1 = F1d / PM * 2.0
W2 = W1 * (1. + Rc) / (1. + nu**2 / qsq)
endif
! For nuclei
IF(IA.gt.2) then
sigt = 0.
sigl = 0.
F1d = 0.
F1p = 0.
! Modifed to use Superscaling from Sick, Donnelly, Maieron,
! nucl-th/0109032
if(IA.eq.2) kf=0.085
if(iA.eq.2) Es=0.0022
! changed 4/09
if(IA.eq.3) kf=0.115
if(iA.eq.3) Es=0.001
! changed 4/09
if(IA.gt.3) kf=0.19
if(iA.gt.3) Es=0.017
if(IA.gt.7) kf=0.228
if(iA.gt.7) Es=0.020
c changed 5/09
if(iA.gt.7) Es=0.0165
if(IA.gt.16) kf=0.230
if(iA.gt.16) Es=0.025
if(IA.gt.25) kf=0.236
if(iA.gt.25) Es=0.018
if(IA.gt.38) kf=0.241
if(iA.gt.38) Es=0.028
if(IA.gt.55) kf=0.241
if(iA.gt.55) Es=0.023
if(IA.gt.60) kf=0.245
if(iA.gt.60) Es=0.028
! changed 5/09
if(iA.gt.55) Es=0.018
! adjust pf to give right width based on kf
pf = 0.5 * kf
! assume this is 2 * pf * qv
DW2DPF = 2. * qv
dw2des = 2. * (nu + PM)
! switched to using 99 bins!
cc do ism = 1,15
cc fyuse = fy(ism)
cc WSQP = WSQ + XX(ISM) * PF * DW2DPF - es * dw2des
do ism = 1,99
fyuse = fyp(ism)/100.
WSQP = WSQ + XXp(ISM) * PF * DW2DPF - es * dw2des
IF(WSQP.GT. 1.159) THEN
call CHRISTY507(Wsqp,Qsq,F1pp,Rc,sigtp,siglp)
call resmodd(wsqp,qsq,xvald0,F1dp)
F1d = F1d + F1dp * Fyuse
F1p = F1p + F1pp * Fyuse
sigt = sigt + sigtp * Fyuse
sigl = sigl + siglp * Fyuse
ENDIF
ENDDO
Rc = 0.
if(sigt .gt. 0.) Rc = sigl / sigt
W1 = (2. * Z * F1d + (A - 2. * Z) * (2. * F1d - F1p)) / PM
W1= W1*(1.0+P(13)*x+P(14)*x**2+P(15)*x**3+P(16)*x**4+P(17)*x**5)
cc if(W .GT. 1.3)
if(W .GT. 0.0)
> W1=W1*(1.0+(P(20)*W+P(21)*W**2)/(1.0+P(22)*QSQ))**2
CALL MEC2009( Z , A , qsq , wsq , F1M )
W1 = W1 + F1M
if(Wsq .gt.0.0 ) Rc = Rc * ( 1.0 + P(6) + P(23)*A )
W2 = W1 * (1. + Rc) / (1. + nu**2 / qsq)
DEBUG=.FALSE.
IF( W1 .LE. 0.0 .AND. DEBUG ) THEN
write(*,*) 'test = ', Z,A,W,QSQ,x,F1M,W1
write(*,*) 'test1 = ',
> (1.0+(P(20)*W+P(21)*W**2)/(1.0+P(22)*QSQ)),
> (1.0+P(13)*x+P(14)*x**2+P(15)*x**3+P(16)*x**4+P(17)*x**5)
ENDIF
ENDIF
A4 = A
x4 = qsq / 2. / pm / nu
emcfac = fitemc(x4, a4, goodfit)
F1 = pm * W1 * emcfac
F2 = nu * W2 * emcfac
RETURN
END
cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
SUBROUTINE MEC2009(z,a,q2,w2,f1)
! fit to low q2 dip region: purefly empirical
! assume contribution is pure transverse
implicit none
real*8 z,a,q2,w2,f1,am/0.9383/,w,nu
integer i
real*8 pb(20)/
> 0.1023E+02, 0.1052E+01, 0.2485E-01, 0.1455E+01,
> 0.5650E+01,-0.2889E+00, 0.4943E-01,-0.8183E-01,
> -0.7495E+00, 0.8426E+00,-0.2829E+01, 0.1607E+01,
> 0.1733E+00, 0.0000E+00, 0.0000E+00, 0.0000E+00,
> 0.0000E+00, 0.0000E+00, 0.0000E+00, 0.0000E+00/
REAL*8 P(0:23)
COMMON/PARCORR/P
real*8 p18
real*8 x, f1corr
f1 = 0.0
if(w2.le.0.0) return
w = sqrt(w2)
nu = (w2 - am**2 + q2) / 2. / am
x = q2 / (2.0 * am * nu )
if(a.lt.2.5) return
p18 = p(18)
! special case for 3He
if(a.gt.2.5 .and. a.lt.3.5) p18 = 70
! special case for 4He
if(a.gt.3.5 .and. a.lt.4.5) p18 = 170.
! new values for C, Al, Cu
if(a.gt.4.5) p18 = 215.
if(a.gt.20.) p18 = 235.
if(a.gt.50.) p18 = 230.
f1corr = P(0)*exp(-((W-P(1))**2)/(P(2)))/
> ((1.0 + MAX( 0.3 , Q2 ) / P(3) ) ** P(4) )*nu**P(5)
> *( 1.0 + P18 * A ** ( 1.0 + P(19) * x ) )
f1 = f1corr
if(f1 .le.1.0E-9 ) f1=0.0
c write(*,*) 'vahe1= ', A, W*W, Q2, f1corr
return
end
SUBROUTINE F1F2QE09(Z, A, qsq, wsq, F1, F2)
c
C Calculates quasielastic A(e,e')X structure functions F1 and F2 PER NUCLEUS
c for A>2 uses superscaling from Sick, Donnelly, Maieron, nucl-th/0109032
c for A=2 uses pre-integrated Paris wave function (see ~bosted/smear.f)
c coded by P. Bosted August to October, 2006
c
c input: Z, A (real*8) Z and A of nucleus (shoud be 2.0D0 for deueron)
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
c
c Note: Deuteron agrees well with Laget (see ~bosted/eg1b/laget.f) for
c a) Q2<1 gev**2 and dsig > 1% of peak: doesnt describe tail at high W
c b) Q2>1 gev**2 on wings of q.e. peak. But, this model is up
c to 50% too big at top of q.e. peak. BUT, F2 DOES agree very
c nicely with Osipenko et al data from CLAS, up to 5 GeV**2
IMPLICIT NONE
REAL*8 P(0:23)
COMMON/PARCORR/P
REAL*8 Z, A, avgN, F1, F2, wsq, qsq
REAL*8 amp/0.93828/, amd/1.8756/
REAL*8 PAULI_SUP1, PAULI_SUP2
REAL*8 GEP, GEN, GMP, GMN, Q, Q3, Q4
REAL*8 RMUP/ 2.792782/ ,RMUN/ -1.913148 /
real*8 pz, nu, dpz, pznom, pzmin
real*8 qv, TAU, W1, W2, FY, dwmin, w2p
real kappa, lam, lamp, taup, squigglef, psi, psip, nuL, nut
real kf, es, GM2bar, GE2bar, W1bar, W2bar, Delta, GL, GT
integer IA, izz, izzmin, izp, izznom, izdif
c Look up tables for deuteron case
real*8 fyd(200)/
> 0.00001,0.00002,0.00003,0.00005,0.00006,0.00009,0.00010,0.00013,
> 0.00015,0.00019,0.00021,0.00026,0.00029,0.00034,0.00038,0.00044,
> 0.00049,0.00057,0.00062,0.00071,0.00078,0.00089,0.00097,0.00109,
> 0.00119,0.00134,0.00146,0.00161,0.00176,0.00195,0.00211,0.00232,
> 0.00252,0.00276,0.00299,0.00326,0.00352,0.00383,0.00412,0.00447,
> 0.00482,0.00521,0.00560,0.00603,0.00648,0.00698,0.00747,0.00803,
> 0.00859,0.00921,0.00985,0.01056,0.01126,0.01205,0.01286,0.01376,
> 0.01467,0.01569,0.01671,0.01793,0.01912,0.02049,0.02196,0.02356,
> 0.02525,0.02723,0.02939,0.03179,0.03453,0.03764,0.04116,0.04533,
> 0.05004,0.05565,0.06232,0.07015,0.07965,0.09093,0.10486,0.12185,
> 0.14268,0.16860,0.20074,0.24129,0.29201,0.35713,0.44012,0.54757,
> 0.68665,0.86965,1.11199,1.43242,1.86532,2.44703,3.22681,4.24972,
> 5.54382,7.04016,8.48123,9.40627,9.40627,8.48123,7.04016,5.54382,
> 4.24972,3.22681,2.44703,1.86532,1.43242,1.11199,0.86965,0.68665,
> 0.54757,0.44012,0.35713,0.29201,0.24129,0.20074,0.16860,0.14268,
> 0.12185,0.10486,0.09093,0.07965,0.07015,0.06232,0.05565,0.05004,
> 0.04533,0.04116,0.03764,0.03453,0.03179,0.02939,0.02723,0.02525,
> 0.02356,0.02196,0.02049,0.01912,0.01793,0.01671,0.01569,0.01467,
> 0.01376,0.01286,0.01205,0.01126,0.01056,0.00985,0.00921,0.00859,
> 0.00803,0.00747,0.00698,0.00648,0.00603,0.00560,0.00521,0.00482,
> 0.00447,0.00412,0.00383,0.00352,0.00326,0.00299,0.00276,0.00252,
> 0.00232,0.00211,0.00195,0.00176,0.00161,0.00146,0.00134,0.00119,
> 0.00109,0.00097,0.00089,0.00078,0.00071,0.00062,0.00057,0.00049,
> 0.00044,0.00038,0.00034,0.00029,0.00026,0.00021,0.00019,0.00015,
> 0.00013,0.00010,0.00009,0.00006,0.00005,0.00003,0.00002,0.00001/
real*8 avp2(200)/
> 1.0,0.98974,0.96975,0.96768,0.94782,0.94450,0.92494,0.92047,
> 0.90090,0.89563,0.87644,0.87018,0.85145,0.84434,0.82593,0.81841,
> 0.80021,0.79212,0.77444,0.76553,0.74866,0.73945,0.72264,0.71343,
> 0.69703,0.68740,0.67149,0.66182,0.64631,0.63630,0.62125,0.61154,
> 0.59671,0.58686,0.57241,0.56283,0.54866,0.53889,0.52528,0.51581,
> 0.50236,0.49291,0.47997,0.47063,0.45803,0.44867,0.43665,0.42744,
> 0.41554,0.40656,0.39511,0.38589,0.37488,0.36611,0.35516,0.34647,
> 0.33571,0.32704,0.31656,0.30783,0.29741,0.28870,0.27820,0.26945,
> 0.25898,0.25010,0.23945,0.23023,0.21943,0.20999,0.19891,0.18911,
> 0.17795,0.16793,0.15669,0.14667,0.13553,0.12569,0.11504,0.10550,
> 0.09557,0.08674,0.07774,0.06974,0.06184,0.05484,0.04802,0.04203,
> 0.03629,0.03129,0.02654,0.02247,0.01867,0.01545,0.01251,0.01015,
> 0.00810,0.00664,0.00541,0.00512,0.00512,0.00541,0.00664,0.00810,
> 0.01015,0.01251,0.01545,0.01867,0.02247,0.02654,0.03129,0.03629,
> 0.04203,0.04802,0.05484,0.06184,0.06974,0.07774,0.08674,0.09557,
> 0.10550,0.11504,0.12569,0.13553,0.14667,0.15669,0.16793,0.17795,
> 0.18911,0.19891,0.20999,0.21943,0.23023,0.23945,0.25010,0.25898,
> 0.26945,0.27820,0.28870,0.29741,0.30783,0.31656,0.32704,0.33571,
> 0.34647,0.35516,0.36611,0.37488,0.38589,0.39511,0.40656,0.41554,
> 0.42744,0.43665,0.44867,0.45803,0.47063,0.47997,0.49291,0.50236,
> 0.51581,0.52528,0.53889,0.54866,0.56283,0.57241,0.58686,0.59671,
> 0.61154,0.62125,0.63630,0.64631,0.66182,0.67149,0.68740,0.69703,
> 0.71343,0.72264,0.73945,0.74866,0.76553,0.77444,0.79212,0.80021,
> 0.81841,0.82593,0.84434,0.85145,0.87018,0.87644,0.89563,0.90090,
> 0.92047,0.92494,0.94450,0.94782,0.96768,0.96975,0.98974,1.0/
c Peter Bosted's correction params
real*8 pb(20)/ 0.1023E+02, 0.1052E+01, 0.2485E-01, 0.1455E+01,
> 0.5650E+01,-0.2889E+00, 0.4943E-01,-0.8183E-01,
> -0.7495E+00, 0.8426E+00,-0.2829E+01, 0.1607E+01,
> 0.1733E+00, 0.0000E+00, 0.0000E+00, 0.0000E+00,
> 0.0000E+00, 0.0000E+00, 0.0000E+00, 0.0000E+00/
real*8 y,R
! return if proton: future change this to allow for
! equivalent W resolution
F1 = 0.
F2 = 0.
IA = int(A)
avgN = A - Z
IF (iA.EQ.1) RETURN
! some kinematic factors. Return if nu or qsq is negative
Nu = (wsq - amp**2 + qsq) / 2. / amp
if(nu .le. 0.0 .or. qsq .lt. 0.) return
TAU = QSQ / 4.0 / amp**2
qv = sqrt(nu**2 + qsq)
! Bosted fit for nucleon form factors Phys. Rev. C 51, p. 409 (1995)
Q = sqrt(QSQ)
Q3 = QSQ * Q
Q4 = QSQ**2
GEP = 1./ (1. + 0.14 * Q + 3.01 * QSQ + 0.02 * Q3 +
> 1.20 * Q4 + 0.32 * Q**5)
GMP = RMUP * GEP
GMN = RMUN / (1.- 1.74 * Q + 9.29 * QSQ - 7.63 * Q3 +
> 4.63 * Q4)
GEN = 1.25 * RMUN * TAU / (1. + 18.3 * TAU) /
> (1. + QSQ / 0.71)**2
! Get kf and Es from superscaling from Sick, Donnelly, Maieron,
c nucl-th/0109032
if(IA.eq.2) kf=0.085
if(iA.eq.2) Es=0.0022
! changed 4/09
if(IA.eq.3) kf=0.115
if(iA.eq.3) Es=0.001
! changed 4/09
if(IA.gt.3) kf=0.19
if(iA.gt.3) Es=0.017
if(IA.gt.7) kf=0.228
if(iA.gt.7) Es=0.020
c changed 5/09
if(iA.gt.7) Es=0.0165
if(IA.gt.16) kf=0.230
if(iA.gt.16) Es=0.025
if(IA.gt.25) kf=0.236
if(iA.gt.25) Es=0.018
if(IA.gt.38) kf=0.241
if(iA.gt.38) Es=0.028
if(IA.gt.55) kf=0.241
if(iA.gt.55) Es=0.023
if(IA.gt.60) kf=0.245
if(iA.gt.60) Es=0.028
! changed 5/09
if(iA.gt.55) Es=0.018
! Pauli suppression model from Tsai RMP 46,816(74) eq.B54
IF((QV .GT. 2.* kf).OR.(iA.EQ.1)) THEN
PAULI_SUP2 =1.0
ELSE
PAULI_SUP2 = 0.75 * (QV / kf) * (1.0 - ((QV / kf)**2)/12.)
ENDIF
PAULI_SUP1 = PAULI_SUP2
! structure functions with off shell factors
kappa = qv / 2. / amp
lam = nu / 2. / amp
lamp = lam - Es / 2. / amp
taup = kappa**2 - lamp**2
squigglef = sqrt(1. + (kf/amp)**2) -1.
! Very close to treshold, could have a problem
if(1.+lamp.le.0.) return
if(taup * (1. + taup).le.0.) return
psi = (lam - tau ) / sqrt(squigglef) /
> sqrt((1.+lam )* tau + kappa * sqrt(tau * (1. + tau)))
psip = (lamp - taup) / sqrt(squigglef) /
> sqrt((1.+lamp)*taup + kappa * sqrt(taup * (1. + taup)))
nuL = (tau / kappa**2)**2
c changed definition of nuT from
c nuT = tau / 2. / kappa**2 + tan(thr/2.)**2
c to this, in order to separate out F1 and F2 (F1 prop. to tan2 term)
nuT = tau / 2. / kappa**2
GM2bar = Pauli_sup1 * (Z * GMP**2 + avgN * GMN**2)
GE2bar = Pauli_sup2 * (Z * GEP**2 + avgN * GEN**2)
W1bar = tau * GM2bar
W2bar = (GE2bar + tau * GM2bar) / (1. + tau)
Delta = squigglef * (1. - psi**2) * (
> sqrt(tau * (1.+tau)) / kappa + squigglef/3. *
> (1. - psi**2) * tau / kappa**2)
GL = kappa**2 / tau * (GE2bar + Delta * W2bar) /
> 2. / kappa / (1. + squigglef * (1. + psi**2) / 2.)
GT = (2. * tau * GM2bar + Delta * W2bar) /
> 2. / kappa / (1. + squigglef * (1. + psi**2) / 2.)
c added to prevent negative xsections:
gt = max(0., gt)
! from Maria Barbaro: see Amaro et al., PRC71,015501(2005).
FY = 1.5576 / (1. + 1.7720**2 * (psip + 0.3014)**2) /
> (1. + exp(-2.4291 * psip)) / kf
! Use PWIA and Paris W.F. for deuteron to get better FY
if(IA.eq.2) then
! value assuming average p2=0.
pz = (qsq - 2. * amp * nu ) / 2. / qv
izz = int((pz + 1.0) / 0.01) + 1
izz = min(200,max(1,izz))
izznom = izz
! ignoring energy term, estimate change in pz to compensate
! for avp2 term
dpz = avp2(izznom) / 2. / qv
izdif = dpz * 150.
dwmin=1.E6
izzmin=0
do izp = izznom, min(200, max(1, izznom + izdif))
pz = -1. + 0.01 * (izp-0.5)
c *** this version gives worse agreement with laget than
c w2p = (amd + nu - sqrt(amp**2 + avp2(izp)))**2 -
c > qv**2 + 2. * qv * pz - avp2(izp)
c this version!
w2p = (amd + nu - amp )**2 -
> qv**2 + 2. * qv * pz - avp2(izp)
c if passed first minimum, quit looking so don't find second one
if(abs(w2p - amp**2).gt.dwmin) goto 11
if(abs(w2p - amp**2).lt.dwmin) then
dwmin = abs(w2p - amp**2)
izzmin = izp
endif
enddo
11 izz = min(199,max(2,izzmin))
! search for minimum in 1/10th bins locally
pznom = -1. + 0.01 * (izz-0.5)
dwmin=1.E6
do izp = 1,19
pz = pznom - 0.01 + 0.001 * izp
c *** this version gives worse agreement with laget than
c w2p = (amd + nu - sqrt(amp**2 + avp2(izz)))**2 -
c > qv**2 + 2. * qv * pz - avp2(izz)
c this version!
w2p = (amd + nu - amp )**2 -
> qv**2 + 2. * qv * pz - avp2(izz)
if(abs(w2p - amp**2).lt.dwmin) then
dwmin = abs(w2p - amp**2)
pzmin = pz
endif
enddo
if(dwmin.ge.1.e6.or.abs(pznom-pzmin).gt.0.01)
> write(6,'(1x,''error in dwmin,pzmin'',3i4,6f7.3)')
> izznom,izzmin,izz,qsq,wsq,w2p,dwmin/1.e6,pzmin,pznom
if(pzmin.lt.pznom) then
fy = fyd(izz) - (fyd(izz-1) - fyd(izz)) *
> (pzmin - pznom) / 0.01
else
fy = fyd(izz) + (fyd(izz+1) - fyd(izz)) *
> (pzmin - pznom) / 0.01
endif
endif
c final results
F2 = nu * FY * (nuL * GL + nuT * GT)
F1 = amp * FY * GT / 2.
if(F1.LT.0.0) F1 = 0.
if(nu.gt.0. .and.f1.gt.0.) then
R = (F2 / nu) / (F1 / amp) * (1. + nu**2 / qsq) - 1.0
else
r = 0.4/qsq
endif
c apply correction factors
if(A.gt.2) then
y = (wsq -amp**2) / qv
c F1 = F1 * (1. + pb(8) + pb(9) * y +
c > pb(10)*y**2 + pb(11)*y**3 + pb(12)*y**4 )
c R = R * (1. + pb(13))
c F2 = nu * F1/amp * (1. + R) / (1. + nu**2/qsq)
cc correction to correction Vahe
if(wsq.gt.0.0) then
F1=F1*(1.0+P(7)+P(8)*y+P(9)*y**2 +P(10)*y**3 +P(11)*y**4)
R = R * ( 1.0 + P(12) )
F2 = nu * F1/amp * (1. + R) / (1. + nu**2/qsq)
if(F1.LT.0.0) F1=0.0
endif
endif
return
end
BLOCK DATA CORRECTION
REAL*8 P(0:23)
COMMON/PARCORR/P
c DATA P/
c c 5.1141e-02, 9.9343e-01, 5.3317e-02, 1.3949e+00,
c c 5.9993e+00, -2.9393e-01, 9.9316e-02, 1.0935e-02,
c c 3.7697e-01, 1.3542e+00, -7.4618e+00, 4.3540e+00,
c c -3.7911e-01, 3.5105e-01, -1.8903e+00, 4.9139e+00,
c c -5.9923e+00, 2.5021e+00, 1.9943e+01, -5.5879e-02,
c c 3.1914e-01, -1.8657e-01, 3.2746e+01, 4.9801e-03/
DATA P/
c 5.1377e-03, 9.8071e-01, 4.6379e-02, 1.6433e+00,
c 6.9826e+00, -2.2655e-01, 1.1095e-01, 2.7945e-02,
c 4.0643e-01, 1.6076e+00, -7.5460e+00, 4.4418e+00,
c -3.7464e-01, 1.0414e-01, -2.6852e-01, 9.6653e-01,
c -1.9055e+00, 9.8965e-01, 2.0613e+02, -4.5536e-02,
c 2.4902e-01, -1.3728e-01, 2.9201e+01, 4.9280e-03/
END
subroutine pind(W2,Q2,F1,R,sigt,sigl)
! Calculate proton with Fermi smearing of a deuteron
implicit none
real*8 q2,w2,F1,R,sigt,sigl,am/0.9383/,nu,qv,F1p,Rp,sigtp,siglp
real*8 amd/1.8756/,w2p,pz
integer ism
real*8 fyd(20)/ 0.4965, 0.4988, 0.4958, 0.5008, 0.5027, 0.5041,
> 0.5029, 0.5034, 0.4993, 0.5147, 0.5140, 0.4975, 0.5007, 0.4992,
> 0.4994, 0.4977, 0.5023, 0.4964, 0.4966, 0.4767/
real*8 avpz(20)/-0.1820,-0.0829,-0.0590,-0.0448,-0.0345,-0.0264,
> -0.0195,-0.0135,-0.0079,-0.0025, 0.0029, 0.0083, 0.0139, 0.0199,
> 0.0268, 0.0349, 0.0453, 0.0598, 0.0844, 0.1853/
real*8 avp2(20)/ 0.0938, 0.0219, 0.0137, 0.0101, 0.0081, 0.0068,
> 0.0060, 0.0054, 0.0051, 0.0049, 0.0050, 0.0051, 0.0055, 0.0060,
> 0.0069, 0.0081, 0.0102, 0.0140, 0.0225, 0.0964/
c Look up tables for deuteron in fine bins for sub threshold
real*8 fydf(200)/
> 0.00001,0.00002,0.00003,0.00005,0.00006,0.00009,0.00010,0.00013,
> 0.00015,0.00019,0.00021,0.00026,0.00029,0.00034,0.00038,0.00044,
> 0.00049,0.00057,0.00062,0.00071,0.00078,0.00089,0.00097,0.00109,
> 0.00119,0.00134,0.00146,0.00161,0.00176,0.00195,0.00211,0.00232,
> 0.00252,0.00276,0.00299,0.00326,0.00352,0.00383,0.00412,0.00447,
> 0.00482,0.00521,0.00560,0.00603,0.00648,0.00698,0.00747,0.00803,
> 0.00859,0.00921,0.00985,0.01056,0.01126,0.01205,0.01286,0.01376,
> 0.01467,0.01569,0.01671,0.01793,0.01912,0.02049,0.02196,0.02356,
> 0.02525,0.02723,0.02939,0.03179,0.03453,0.03764,0.04116,0.04533,
> 0.05004,0.05565,0.06232,0.07015,0.07965,0.09093,0.10486,0.12185,
> 0.14268,0.16860,0.20074,0.24129,0.29201,0.35713,0.44012,0.54757,
> 0.68665,0.86965,1.11199,1.43242,1.86532,2.44703,3.22681,4.24972,
> 5.54382,7.04016,8.48123,9.40627,9.40627,8.48123,7.04016,5.54382,
> 4.24972,3.22681,2.44703,1.86532,1.43242,1.11199,0.86965,0.68665,
> 0.54757,0.44012,0.35713,0.29201,0.24129,0.20074,0.16860,0.14268,
> 0.12185,0.10486,0.09093,0.07965,0.07015,0.06232,0.05565,0.05004,
> 0.04533,0.04116,0.03764,0.03453,0.03179,0.02939,0.02723,0.02525,
> 0.02356,0.02196,0.02049,0.01912,0.01793,0.01671,0.01569,0.01467,
> 0.01376,0.01286,0.01205,0.01126,0.01056,0.00985,0.00921,0.00859,
> 0.00803,0.00747,0.00698,0.00648,0.00603,0.00560,0.00521,0.00482,
> 0.00447,0.00412,0.00383,0.00352,0.00326,0.00299,0.00276,0.00252,
> 0.00232,0.00211,0.00195,0.00176,0.00161,0.00146,0.00134,0.00119,
> 0.00109,0.00097,0.00089,0.00078,0.00071,0.00062,0.00057,0.00049,
> 0.00044,0.00038,0.00034,0.00029,0.00026,0.00021,0.00019,0.00015,
> 0.00013,0.00010,0.00009,0.00006,0.00005,0.00003,0.00002,0.00001/
real*8 avp2f(200)/
> 1.0,0.98974,0.96975,0.96768,0.94782,0.94450,0.92494,0.92047,
> 0.90090,0.89563,0.87644,0.87018,0.85145,0.84434,0.82593,0.81841,
> 0.80021,0.79212,0.77444,0.76553,0.74866,0.73945,0.72264,0.71343,
> 0.69703,0.68740,0.67149,0.66182,0.64631,0.63630,0.62125,0.61154,
> 0.59671,0.58686,0.57241,0.56283,0.54866,0.53889,0.52528,0.51581,
> 0.50236,0.49291,0.47997,0.47063,0.45803,0.44867,0.43665,0.42744,
> 0.41554,0.40656,0.39511,0.38589,0.37488,0.36611,0.35516,0.34647,
> 0.33571,0.32704,0.31656,0.30783,0.29741,0.28870,0.27820,0.26945,
> 0.25898,0.25010,0.23945,0.23023,0.21943,0.20999,0.19891,0.18911,
> 0.17795,0.16793,0.15669,0.14667,0.13553,0.12569,0.11504,0.10550,
> 0.09557,0.08674,0.07774,0.06974,0.06184,0.05484,0.04802,0.04203,
> 0.03629,0.03129,0.02654,0.02247,0.01867,0.01545,0.01251,0.01015,
> 0.00810,0.00664,0.00541,0.00512,0.00512,0.00541,0.00664,0.00810,
> 0.01015,0.01251,0.01545,0.01867,0.02247,0.02654,0.03129,0.03629,
> 0.04203,0.04802,0.05484,0.06184,0.06974,0.07774,0.08674,0.09557,
> 0.10550,0.11504,0.12569,0.13553,0.14667,0.15669,0.16793,0.17795,
> 0.18911,0.19891,0.20999,0.21943,0.23023,0.23945,0.25010,0.25898,
> 0.26945,0.27820,0.28870,0.29741,0.30783,0.31656,0.32704,0.33571,
> 0.34647,0.35516,0.36611,0.37488,0.38589,0.39511,0.40656,0.41554,
> 0.42744,0.43665,0.44867,0.45803,0.47063,0.47997,0.49291,0.50236,
> 0.51581,0.52528,0.53889,0.54866,0.56283,0.57241,0.58686,0.59671,
> 0.61154,0.62125,0.63630,0.64631,0.66182,0.67149,0.68740,0.69703,
> 0.71343,0.72264,0.73945,0.74866,0.76553,0.77444,0.79212,0.80021,
> 0.81841,0.82593,0.84434,0.85145,0.87018,0.87644,0.89563,0.90090,
> 0.92047,0.92494,0.94450,0.94782,0.96768,0.96975,0.98974,1.0/
nu = (w2 - am**2 + q2) / 2. / am
qv = sqrt(nu**2 + q2)
F1=0.
R=0.
sigt=0.
sigl=0.
! Do fast 20 bins if abvoe threshold
if(w2.gt.1.30) then
do ism = 1,20
w2p = (amd + nu - sqrt(am**2 + avp2(ism)))**2 -
> qv**2 + 2. * qv * avpz(ism) - avp2(ism)
if(w2p.gt.1.155) then
call CHRISTY507(W2p,Q2,F1p,Rp,sigtp,siglp)
sigt = sigt + sigtp * fyd(ism) / 10.
sigl = sigl + siglp * fyd(ism) / 10.
F1 = F1 + F1p * fyd(ism) / 10.
endif
enddo
else
do ism = 1,200
pz = -1. + 0.01 * (ism-0.5)
c Need avp2f term to get right behavior x>1!
w2p = (amd + nu - sqrt(am**2 + avp2f(ism)))**2 -
> qv**2 + 2. * qv * pz - avp2f(ism)
if(w2p.gt.1.155) then
call CHRISTY507(W2p,Q2,F1p,Rp,sigtp,siglp)
sigt = sigt + sigtp * fydf(ism) / 100.
sigl = sigl + siglp * fydf(ism) / 100.
F1 = F1 + F1p * fydf(ism) / 100.
endif
enddo
endif
if(sigt.ne.0.) R = sigl / sigt
return
end
subroutine resd(q2,w2,xval,F1)
! Calculate dueteron F1 by Fermi smearing of proton plus neutron
! Add MEC term (not smeared)
c P.E. Bosted and M.E. Christy, ``Empirical Fit to Inelastic
c Electron-Deuteron and Electron-Neutron
c Resonance Region Transverse Cross Sections,
c (arXiv:0711.0159). Submitted to Phys. Rev. C.
implicit none
real*8 q2,w2,xval(50),F1,am/0.9383/,nu,qv,dw2dpf,w2p,sigp,f1sv
real*8 sigtst,amd/1.8756/, pz, f1m,f2m
integer ism,i
real*8 fyd(20)/ 0.4965, 0.4988, 0.4958, 0.5008, 0.5027, 0.5041,
> 0.5029, 0.5034, 0.4993, 0.5147, 0.5140, 0.4975, 0.5007, 0.4992,
> 0.4994, 0.4977, 0.5023, 0.4964, 0.4966, 0.4767/
real*8 avpz(20)/-0.1820,-0.0829,-0.0590,-0.0448,-0.0345,-0.0264,
> -0.0195,-0.0135,-0.0079,-0.0025, 0.0029, 0.0083, 0.0139, 0.0199,
> 0.0268, 0.0349, 0.0453, 0.0598, 0.0844, 0.1853/
real*8 avp2(20)/ 0.0938, 0.0219, 0.0137, 0.0101, 0.0081, 0.0068,
> 0.0060, 0.0054, 0.0051, 0.0049, 0.0050, 0.0051, 0.0055, 0.0060,
> 0.0069, 0.0081, 0.0102, 0.0140, 0.0225, 0.0964/
c Look up tables for deuteron in fine bins for sub threshold
real*8 fydf(200)/
> 0.00001,0.00002,0.00003,0.00005,0.00006,0.00009,0.00010,0.00013,
> 0.00015,0.00019,0.00021,0.00026,0.00029,0.00034,0.00038,0.00044,
> 0.00049,0.00057,0.00062,0.00071,0.00078,0.00089,0.00097,0.00109,
> 0.00119,0.00134,0.00146,0.00161,0.00176,0.00195,0.00211,0.00232,
> 0.00252,0.00276,0.00299,0.00326,0.00352,0.00383,0.00412,0.00447,
> 0.00482,0.00521,0.00560,0.00603,0.00648,0.00698,0.00747,0.00803,
> 0.00859,0.00921,0.00985,0.01056,0.01126,0.01205,0.01286,0.01376,
> 0.01467,0.01569,0.01671,0.01793,0.01912,0.02049,0.02196,0.02356,
> 0.02525,0.02723,0.02939,0.03179,0.03453,0.03764,0.04116,0.04533,
> 0.05004,0.05565,0.06232,0.07015,0.07965,0.09093,0.10486,0.12185,
> 0.14268,0.16860,0.20074,0.24129,0.29201,0.35713,0.44012,0.54757,
> 0.68665,0.86965,1.11199,1.43242,1.86532,2.44703,3.22681,4.24972,
> 5.54382,7.04016,8.48123,9.40627,9.40627,8.48123,7.04016,5.54382,
> 4.24972,3.22681,2.44703,1.86532,1.43242,1.11199,0.86965,0.68665,
> 0.54757,0.44012,0.35713,0.29201,0.24129,0.20074,0.16860,0.14268,
> 0.12185,0.10486,0.09093,0.07965,0.07015,0.06232,0.05565,0.05004,
> 0.04533,0.04116,0.03764,0.03453,0.03179,0.02939,0.02723,0.02525,
> 0.02356,0.02196,0.02049,0.01912,0.01793,0.01671,0.01569,0.01467,
> 0.01376,0.01286,0.01205,0.01126,0.01056,0.00985,0.00921,0.00859,
> 0.00803,0.00747,0.00698,0.00648,0.00603,0.00560,0.00521,0.00482,
> 0.00447,0.00412,0.00383,0.00352,0.00326,0.00299,0.00276,0.00252,
> 0.00232,0.00211,0.00195,0.00176,0.00161,0.00146,0.00134,0.00119,
> 0.00109,0.00097,0.00089,0.00078,0.00071,0.00062,0.00057,0.00049,
> 0.00044,0.00038,0.00034,0.00029,0.00026,0.00021,0.00019,0.00015,
> 0.00013,0.00010,0.00009,0.00006,0.00005,0.00003,0.00002,0.00001/
real*8 avp2f(200)/
> 1.0,0.98974,0.96975,0.96768,0.94782,0.94450,0.92494,0.92047,
> 0.90090,0.89563,0.87644,0.87018,0.85145,0.84434,0.82593,0.81841,
> 0.80021,0.79212,0.77444,0.76553,0.74866,0.73945,0.72264,0.71343,
> 0.69703,0.68740,0.67149,0.66182,0.64631,0.63630,0.62125,0.61154,
> 0.59671,0.58686,0.57241,0.56283,0.54866,0.53889,0.52528,0.51581,
> 0.50236,0.49291,0.47997,0.47063,0.45803,0.44867,0.43665,0.42744,
> 0.41554,0.40656,0.39511,0.38589,0.37488,0.36611,0.35516,0.34647,
> 0.33571,0.32704,0.31656,0.30783,0.29741,0.28870,0.27820,0.26945,
> 0.25898,0.25010,0.23945,0.23023,0.21943,0.20999,0.19891,0.18911,
> 0.17795,0.16793,0.15669,0.14667,0.13553,0.12569,0.11504,0.10550,
> 0.09557,0.08674,0.07774,0.06974,0.06184,0.05484,0.04802,0.04203,
> 0.03629,0.03129,0.02654,0.02247,0.01867,0.01545,0.01251,0.01015,
> 0.00810,0.00664,0.00541,0.00512,0.00512,0.00541,0.00664,0.00810,
> 0.01015,0.01251,0.01545,0.01867,0.02247,0.02654,0.03129,0.03629,
> 0.04203,0.04802,0.05484,0.06184,0.06974,0.07774,0.08674,0.09557,
> 0.10550,0.11504,0.12569,0.13553,0.14667,0.15669,0.16793,0.17795,
> 0.18911,0.19891,0.20999,0.21943,0.23023,0.23945,0.25010,0.25898,
> 0.26945,0.27820,0.28870,0.29741,0.30783,0.31656,0.32704,0.33571,
> 0.34647,0.35516,0.36611,0.37488,0.38589,0.39511,0.40656,0.41554,
> 0.42744,0.43665,0.44867,0.45803,0.47063,0.47997,0.49291,0.50236,
> 0.51581,0.52528,0.53889,0.54866,0.56283,0.57241,0.58686,0.59671,
> 0.61154,0.62125,0.63630,0.64631,0.66182,0.67149,0.68740,0.69703,
> 0.71343,0.72264,0.73945,0.74866,0.76553,0.77444,0.79212,0.80021,
> 0.81841,0.82593,0.84434,0.85145,0.87018,0.87644,0.89563,0.90090,
> 0.92047,0.92494,0.94450,0.94782,0.96768,0.96975,0.98974,1.0/
logical usemec/.true./
nu = (w2 - am**2 + q2) / 2. / am
qv = sqrt(nu**2 + q2)
F1 = 0.
! Do fast 20 bins if abvoe threshold
if(w2.gt.1.30) then
do ism = 1,20
w2p = (amd + nu - sqrt(am**2 + avp2(ism)))**2 -
> qv**2 + 2. * qv * avpz(ism) - avp2(ism)
if(w2p.gt.1.155) then
call resmodd(w2p,q2,xval,sigp)
F1 = F1 + sigp * fyd(ism) / 10.
endif
enddo
else
f1 = 0.
do ism = 1,200
pz = -1. + 0.01 * (ism-0.5)
! Need avp2f term to get right behavior x>1!
w2p = (amd + nu - sqrt(am**2 + avp2f(ism)))**2 -
> qv**2 + 2. * qv * pz - avp2f(ism)
if(w2p.gt.1.155) then
call resmodd(w2p,q2,xval,sigp)
F1 = F1 + sigp * fydf(ism) / 100.
endif
enddo
endif
! add MEC term
! took out again, then put back in again
call mec(1.D0,2.D0,q2,w2,f1m,f2m,xval)
if(usemec) f1 = f1 + f1m
return
end
subroutine resder(q2,w2,xval,F1er)
! Get error on F1 from deuteron fit
implicit none
integer i,j,icall
real*8 q2,w2,xval(50),xvalp(50),deriv(50),emat(50,50),F1,F1er
real*8 epsilon,F1p
integer ifix(50)/1, 2, 3, 4, 5, 6, 7, 0, 0, 0,
> 0, 0, 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, 0, 0,42,0/
logical first/.true./
icall = icall + 1
if(first) then
open(unit=9,file='F1F207D2emat.dat')
do i=1,50
read(9,'(10E12.4)') (emat(i,j),j=1,50)
enddo
close(unit=9)
first = .false.
endif
call resd(q2,w2,xval,F1)
epsilon = 1.E-7
do i=1,50
do j=1,50
xvalp(j) = xval(j)
if(i.eq.j) xvalp(j) = xvalp(j) + epsilon
enddo
call resd(q2,w2,xvalp,F1p)
deriv(i) = (F1p - F1) / epsilon
if(icall.eq.10) write(6,'(i3,3e12.4)') i,f1,f1p,deriv(i)
enddo
F1er=0.
do i=1,50
do j=1,50
if(ifix(i).ne.0 .and.ifix(j).ne.0) then
F1er = F1er + emat(ifix(i),ifix(j)) * deriv(i) * deriv(j)
endif
enddo
enddo
! add overall Hall C norm. err
F1er = sqrt(F1er + (0.02*F1)**2)
return
end
! returns F1 for average of free proton and neutron
! for given W2, Q2
SUBROUTINE RESMODD(w2,q2,xval,sig)
IMPLICIT NONE
REAL*8 W,w2,q2,mp,mp2,mpi2,xb,xth(4),sig,xval(50),mass(7),width(7)
REAL*8 height(7),sig_del,sig_21,sig_22,sig_31,sig_32,rescoef(6,4)
REAL*8 nr_coef(3,4),sigr(5000,7),wdif,wdif2,wr,sig_nr,sig_4,q2temp
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,2),spin(7),ang(7)
REAL*8 pgam(7),pwid(7,2),x0(7),dip,dip2
REAL*8 sig_res,sig_4L,xpr,alpha,pi,F1
INTEGER i,j,l,lmax,num,iw
real*8 noverp,fnp_nmc,x,a,b,sig_mec,brp(7,3)
real*8 xval0(12)/
c new 5/07 values. Values 1 and 7 will be overridden below.
& 0.12298E+01,0.15304E+01,0.15057E+01,0.16980E+01,0.16650E+01,
& 0.14333E+01,0.13573E+00,0.22000E+00,0.82956E-01,0.95782E-01,
& 0.10936E+00,0.37944E+00/
real*8 xvalold(50),w2sv,q2sv,sigrsv(7),md,w2p,wp,wdifp,xprp,nu
logical first/.true./
common/tst2/sigrsv,sig_nr,sig_mec
real br2(7),br3(7)
save
sig = 0.
if(w2.lt.1.07327**2 .or. w2.gt.25 .or.
> q2.lt.0.0 .or. q2.gt.11.0) then
write(6,'(1x,''error, q2 or w2 out of range'',2f8.3)') w2,q2
return
endif
c do this if fitting masses or widths, else set first true in above
first=.false.
if(xvalold(1).ne.xval(1) .or.
> xvalold(7).ne.xval(7)) first=.true.
xvalold(1)=xval(1)
xvalold(7)=xval(7)
if(first) then
mp = 0.9382727
mpi = 0.135
mpi2 = mpi*mpi
meta = 0.547
mp2 = mp*mp
pi = 3.141593
alpha = 1./137.036
! branching ratios
br(1,1) = 1.0
br(2,1) = 0.5
br(3,1) = 0.65
br(4,1) = 0.65
br(5,1) = 0.4
br(6,1) = 0.65
br(7,1) = 0.6
! angular momenta
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. !!! ? 4th resonance region
! x0 parameter
do i=1,7
c 2006
c x0(i) = 0.165