00001 #ifndef __CWAVEFUNCTION_WFCOMMON_CC__
00002 #define __CWAVEFUNCTION_WFCOMMON_CC__
00003 #include "wavefunction.h"
00004
00005 #include <cstdio>
00006 #include <cstdlib>
00007 #include <cmath>
00008 #include <complex>
00009 #include <cstring>
00010
00011 using namespace std;
00012
00013 CWaveFunction::CWaveFunction(){
00014 generic=0;
00015 ci=complex<double>(0.0,1.0);
00016 MPI=139.58;
00017 MKAON=493.677;
00018 MPROTON=938.271;
00019 MNEUTRON=939.565;
00020 MLAMBDA=1115.7;
00021 }
00022
00023 CWaveFunction::~CWaveFunction(){
00024 int ichannel,iq,ell;
00025 for(ichannel=0;ichannel<nchannels;ichannel++){
00026 delete Wepsilon[ichannel];
00027 delete delta[ichannel];
00028 delete ddeltadq[ichannel];
00029 }
00030 delete [] Wepsilon;
00031 delete [] delta;
00032 delete [] ddeltadq;
00033 delete [] qarray;
00034 delete [] eta;
00035 delete [] channelweight;
00036
00037 for(iq=0;iq<nqmax;iq++) delete (planewave[iq]);
00038 delete [] planewave;
00039 for(ell=0;ell<=ellmax;ell++){
00040 for(iq=0;iq<nqmax;iq++){
00041 if(partwave[ell][iq]!=NULL) delete (partwave[ell][iq]);
00042 }
00043 delete [] partwave[ichannel];
00044 }
00045 delete [] partwave;
00046 }
00047
00048 bool CWaveFunction::GetIDENTICAL(){
00049 return IDENTICAL;
00050 }
00051
00052 int CWaveFunction::GetNQMAX(){
00053 return nqmax;
00054 }
00055
00056 int CWaveFunction::GetNCHANNELS(){
00057 return nchannels;
00058 }
00059
00060 double CWaveFunction::GetDELTA(int ichannel,int iq){
00061 return delta[ichannel][iq];
00062 }
00063
00064 double CWaveFunction::GetDELQ(){
00065 return delq;
00066 }
00067
00068 double CWaveFunction::GetQ(int iq){
00069 return qarray[iq];
00070 }
00071
00072 void CWaveFunction::ParsInit(string parsfilename){
00073 parameterMap parameters;
00074 FILE *qarrayfile;
00075 char qarrayfilename[120];
00076 string stemp;
00077 int iq;
00078 bool filetest=0;
00079
00080 parameter::ReadParsFromFile(parameters,parsfilename);
00081
00082 stemp=parameter::getS(parameters,"QARRAYFILENAME","no_qarray_file");
00083 strcpy(qarrayfilename,stemp.c_str());
00084 printf("qarrayfilename set to %s\n",qarrayfilename);
00085
00086 delq=parameter::getD(parameters,"DELQ",-999);
00087 nqmax=parameter::getI(parameters,"NQMAX",-999);
00088 epsilon=parameter::getD(parameters,"EPSILON",-999);
00089 COULOMB=parameter::getB(parameters,"COULOMB",-999);
00090 STRONG=parameter::getB(parameters,"STRONG",-999);
00091 IDENTICAL=parameter::getB(parameters,"IDENTICAL",0);
00092
00093
00094 if(delq<0) filetest=1;
00095 if(filetest==1){
00096 parameter::set(parameters,"delq",-1.0);
00097 printf("will read qarray from %s\n",qarrayfilename);
00098 qarrayfile=fopen(qarrayfilename,"r");
00099 fscanf(qarrayfile,"%d",&nqmax);
00100 parameter::set(parameters,"NQMAX",nqmax);
00101 qarray=new double[nqmax];
00102 for(iq=0;iq<nqmax;iq++){
00103 fscanf(qarrayfile,"%lf",&qarray[iq]);
00104 }
00105 fclose(qarrayfile);
00106 }
00107 else{
00108 delq=parameter::getD(parameters,"DELQ",-999);
00109 qarray=new double[nqmax];
00110 for(iq=0;iq<nqmax;iq++){
00111 qarray[iq]=(iq+0.5)*delq;
00112 }
00113 }
00114 printf("__ WFPARAMETERS __\n");
00115 printf("delq set to %g\n",delq);
00116 printf("nqmax set to %d\n",nqmax);
00117 printf("epsilon set to %g\n",epsilon);
00118 printf("STRONG set to %d\n",STRONG);
00119 printf("COULOMB set to %d\n",COULOMB);
00120 printf("IDENTICAL set to %d\n",IDENTICAL);
00121 printf("____________________\n");
00122 }
00123
00124 void CWaveFunction::InitArrays(){
00125 int iq,l,ichannel;
00126 eta=new double[nqmax];
00127 planewave=new CPlaneWave*[nqmax];
00128 partwave=new CPartWave **[ellmax+1];
00129 for(l=0;l<=ellmax;l++){
00130 partwave[l]=new CPartWave *[nqmax];
00131 for(iq=0;iq<nqmax;iq++) partwave[l][iq]=NULL;
00132 }
00133
00134 delta=new double *[nchannels];
00135 ddeltadq=new double *[nchannels];
00136 Wepsilon=new double *[nchannels];
00137 for(ichannel=0;ichannel<nchannels;ichannel++){
00138 delta[ichannel]=new double[nqmax];
00139 ddeltadq[ichannel]=new double[nqmax];
00140 Wepsilon[ichannel]=new double[nqmax];
00141 }
00142
00143 ell=new int[nchannels];
00144 channelweight=new double[nchannels];
00145
00146 }
00147
00148 void CWaveFunction::InitWaves(){
00149 int iq,ichannel;
00150 double q,e1,e2,e;
00151
00152 for(iq=0;iq<nqmax;iq++){
00153 q=qarray[iq];
00154 e1=sqrt(m1*m1+q*q);
00155 e2=sqrt(m2*m2+q*q);
00156 e=e1+e2;
00157
00158 eta[iq]=double(q1q2)*(pow(e,4)-pow(m1*m1-m2*m2,2))/(4.0*e*e*e*137.036*q);
00159
00160 planewave[iq]=new CPlaneWave(eta[iq],q1q2,q);
00161 }
00162 for(ichannel=0;ichannel<nchannels;ichannel++){
00163 for(iq=0;iq<nqmax;iq++){
00164 q=qarray[iq];
00165 if(partwave[ell[ichannel]][iq]==NULL){
00166 partwave[ell[ichannel]][iq]
00167 =new CPartWave(eta[iq],q1q2,q,ell[ichannel],epsilon);
00168 }
00169 }
00170 }
00171 }
00172
00173 double CWaveFunction::GetPsiSquared(double q,double r,double ctheta){
00174 int iq,iqlow,iqhigh;
00175 double wlow,whigh,interpolate,qscaled,rscaled;
00176 if(r<1000.0){
00177 if(generic==1 && q1q2!=0){
00178 qscaled=q*(muscale/mu)*q1q2scale/double(q1q2);
00179 rscaled=q*r/qscaled;
00180 }
00181 else{
00182 qscaled=q;
00183 rscaled=r;
00184 }
00185
00186 if(delq<0){
00187 iq=0;
00188 while(qscaled>qarray[iq]+1.0E-5 && iq<nqmax){
00189 iq+=1;
00190 }
00191 iqlow=iq-1;
00192 iqhigh=iq;
00193 }
00194 else{
00195 iqhigh=lrint(qscaled/delq);
00196 if(iqhigh==0) iqhigh=1;
00197 iqlow=iqhigh-1;
00198 }
00199 if(iqhigh==nqmax && nqmax>1 && (qscaled-qarray[nqmax-1])<0.5*(qarray[nqmax-1]-qarray[nqmax-2])){
00200 iqhigh=nqmax-1;
00201 iqlow=iqhigh-1;
00202 }
00203
00204 if(iqhigh==0) return CalcPsiSquared(0,rscaled,ctheta);
00205 else if(iqhigh>=nqmax && qscaled-qarray[nqmax-1]>1.0E-5) return 1.0;
00206 else{
00207 wlow=(qarray[iqhigh]-qscaled)/(qarray[iqhigh]-qarray[iqlow]);
00208 whigh=1.0-wlow;
00209 if(fabs(wlow)<1.0E-5) interpolate=CalcPsiSquared(iqhigh,rscaled,ctheta);
00210 else if(fabs(whigh)<1.0E-5)
00211 interpolate=CalcPsiSquared(iqlow,rscaled,ctheta);
00212 else{
00213 interpolate=wlow*CalcPsiSquared(iqlow,rscaled,ctheta)
00214 +whigh*CalcPsiSquared(iqhigh,rscaled,ctheta);
00215 }
00216 if(rscaled>1.0 && qscaled>qarray[0] && qscaled<qarray[nqmax-1] && interpolate<-0.01){
00217 printf("interpolate=%g, qscaled=%g, qlow=%g, qhigh=%g, rscaled=%g\n",
00218 interpolate,qscaled,qarray[iqlow],qarray[iqhigh],rscaled);
00219 printf("wlow=%g, whigh=%g\n",wlow,whigh);
00220 printf("iqlow=%d, %g\n",iqlow,CalcPsiSquared(iqlow,rscaled,ctheta));
00221 printf("iqhigh=%d, %g\n",iqhigh,CalcPsiSquared(iqhigh,rscaled,ctheta));
00222
00223 }
00224 return interpolate;
00225 }
00226 }
00227 else return 1.0;
00228 }
00229
00230 double CWaveFunction::GetPsiSquared(double *pa,double *xa,double *pb,double *xb){
00231 double q,r,ctheta;
00232 getqrctheta(pa,xa,pb,xb,&q,&r,&ctheta);
00233 return GetPsiSquared(q,r,ctheta);
00234 }
00235
00236 double CWaveFunction::CalcPsiSquared(int iq,double r,double ctheta){
00237 return 1.0;
00238 }
00239
00240 void CWaveFunction::getqrctheta(double *p1,double *r1,double *p2,double *r2,double *q,double *r,double *ctheta){
00241 int alpha;
00242 const double g[4]={1.0,-1.0,-1.0,-1.0};
00243 double n[4],qvec[4],rvec[4],nnorm,ndotq,ndotr;
00244
00245 nnorm=0.0;
00246 ndotq=0.0;
00247 ndotr=0.0;
00248 for(alpha=0;alpha<4;alpha++){
00249 n[alpha]=p1[alpha]+p2[alpha];
00250 qvec[alpha]=0.5*(p1[alpha]-p2[alpha]);
00251 rvec[alpha]=r1[alpha]-r2[alpha];
00252 nnorm+=g[alpha]*n[alpha]*n[alpha];
00253 ndotq+=n[alpha]*qvec[alpha]*g[alpha];
00254 ndotr+=n[alpha]*rvec[alpha]*g[alpha];
00255 }
00256 nnorm=sqrt(nnorm);
00257 ndotq=ndotq/nnorm;
00258 ndotr=ndotr/nnorm;
00259
00260 *ctheta=0.0;
00261 *r=0.0;
00262 *q=0.0;
00263 for(alpha=0;alpha<4;alpha++){
00264 n[alpha]=n[alpha]/nnorm;
00265 rvec[alpha]=rvec[alpha]-ndotr*n[alpha];
00266 qvec[alpha]=qvec[alpha]-ndotq*n[alpha];
00267 *r-=g[alpha]*rvec[alpha]*rvec[alpha];
00268 *q-=g[alpha]*qvec[alpha]*qvec[alpha];
00269 *ctheta-=g[alpha]*rvec[alpha]*qvec[alpha];
00270 }
00271 if(*r<0.0 || *q<0.0 || fabs(*ctheta)>sqrt(*q**r)){
00272 printf("Disaster, r^2=%g, q^2=%g, ctheta=%g\n",*r,*q,*ctheta/sqrt(*r**q));
00273 exit(1);
00274 }
00275 *r=sqrt(*r);
00276 *q=sqrt(*q);
00277 *ctheta=*ctheta/(*r**q);
00278 }
00279
00280 void CWaveFunction::PrintPhaseShifts(){
00281 int iq,ichannel;
00282 printf("-------- PHASE SHIFTS --------\n");
00283 printf("q(MeV/c)");
00284 for(ichannel=0;ichannel<nchannels;ichannel++) printf(" l=%d ",ell[ichannel]);
00285 printf("\n");
00286 for(iq=0;iq<nqmax;iq++){
00287 printf("%7.2f ",GetQ(iq));
00288 for(ichannel=0;ichannel<nchannels;ichannel++) printf("% 10.3f",(180.0/PI)*delta[ichannel][iq]);
00289 printf("\n");
00290 }
00291 }
00292
00293 void CWaveFunction::PrintCdelta(double Rx,double Ry,double Rz){
00294 double q,clocal;
00295 int ichannel,iq;
00296 printf("! Qinv C(Q)_estimated ~ ddelta/dq\n");
00297 for(iq=0;iq<nqmax;iq++){
00298 q=qarray[iq];
00299 clocal=1.0;
00300 for(ichannel=0;ichannel<nchannels;ichannel++){
00301 clocal+=(channelweight[ichannel]*(2.0*PI)*pow(HBARC,3)
00302 /(q*q*Rx*Ry*Rz*pow(4.0*PI,1.5)))
00303 *ddeltadq[ichannel][iq];
00304 }
00305 printf("%6.2f %8.4f %g\n",q,clocal,4.0*q*q*(clocal-1.0));
00306 }
00307 printf("_________________________________\n");
00308 }
00309
00310 double CWaveFunction::RelativisticCorrection(double r,int iq){
00311 if(q1q2==0) return 1.0;
00312 else{
00313 const double ALPHA=1.0/137.036;
00314 double q,E,E1,E2,dmudE;
00315 q=GetQ(iq);
00316 E1=sqrt(m1*m1+q*q);
00317 E2=sqrt(m2*m2+q*q);
00318 E=E1+E2;
00319 dmudE=0.25*(1.0-(3.0/pow(E,4))*pow(m1*m1-m2*m2,2));
00320 return 1.0-dmudE*(E/(E1*E2))*q1q2*HBARC*ALPHA/r;
00321 }
00322 }
00323
00324 void CWaveFunction::EffectiveRange(int ichannel,double scattlength,double Reff){
00325 int iq;
00326 double q,tandel,a;
00327 a=scattlength;
00328 for(iq=0;iq<nqmax;iq++){
00329 q=qarray[iq];
00330 tandel=1.0/( (-HBARC/(q*a))+0.5*q*Reff/HBARC);
00331 delta[ichannel][iq]=atan(tandel);
00332 ddeltadq[ichannel][iq]=(tandel*tandel/(1.0+tandel*tandel))
00333 *((-HBARC/(q*q*a))-0.5*Reff/HBARC);
00334
00335
00336 }
00337 }
00338 double CWaveFunction::GetIW(int ell,double epsilon,double q,int q1q2,double eta,double delta){
00339 complex<double> psi0,psi,psi2,psiminus0,psiminus,psiminus2;
00340 complex<double> psiplus0,psiplus,psiplus2;
00341 complex<double> ddeta_psi,ddeta_psiminus,ddeta_psiplus;
00342 complex<double> I,e2idelta;
00343 double x,deleta,a2,root;
00344
00345 if(q1q2!=0){
00346 deleta=0.002*eta;
00347 x=q*epsilon/HBARC;
00348 e2idelta=Misc::ceiphi(2.0*delta);
00349 if(ell==0){
00350 a2=1.0+eta*eta;
00351 root=sqrt(a2);
00352
00353 psi0=CoulWave::CWoutgoing(ell,x,eta-0.5*deleta);
00354 psi0=psi0*e2idelta+conj(psi0);
00355
00356 psiplus0=CoulWave::CWoutgoing(ell+1,x,eta-0.5*deleta);
00357 psiplus0=psiplus0*e2idelta+conj(psiplus0);
00358
00359 psi=CoulWave::CWoutgoing(ell,x,eta);
00360 psi=psi*e2idelta+conj(psi);
00361
00362 psiplus=CoulWave::CWoutgoing(ell+1,x,eta);
00363 psiplus=psiplus*e2idelta+conj(psiplus);
00364
00365 psi2=CoulWave::CWoutgoing(ell,x,eta+0.5*deleta);
00366 psi2=psi2*e2idelta+conj(psi2);
00367
00368 psiplus2=CoulWave::CWoutgoing(ell+1,x,eta+0.5*deleta);
00369 psiplus2=psiplus2*e2idelta+conj(psiplus2);
00370
00371 ddeta_psi=(psi2-psi0)/deleta;
00372 ddeta_psiplus=(psiplus2-psiplus0)/deleta;
00373
00374 I=(conj(psi)*psi+conj(psiplus)*psiplus)*x*a2;
00375 I-=conj(psi)*psiplus*(a2*(1.0+2.0*eta*x)+eta*eta)/root;
00376 I+=(conj(psiplus)*ddeta_psi-conj(psi)*ddeta_psiplus)*eta*root;
00377 }
00378 else{
00379 a2=double(ell*ell)+eta*eta;
00380 root=sqrt(a2);
00381
00382 psi0=CoulWave::CWoutgoing(ell,x,eta-0.5*deleta);
00383 psi0=psi0*e2idelta+conj(psi0);
00384
00385 psiminus0=CoulWave::CWoutgoing(ell-1,x,eta-0.5*deleta);
00386 psiminus0=psiminus0*e2idelta+conj(psiminus0);
00387
00388 psi=CoulWave::CWoutgoing(ell,x,eta);
00389 psi=psi*e2idelta+conj(psi);
00390
00391 psiminus=CoulWave::CWoutgoing(ell-1,x,eta);
00392 psiminus=psiminus*e2idelta+conj(psiminus);
00393
00394 psi2=CoulWave::CWoutgoing(ell,x,eta+0.5*deleta);
00395 psi2=psi2*e2idelta+conj(psi2);
00396
00397 psiminus2=CoulWave::CWoutgoing(ell-1,x,eta+0.5*deleta);
00398 psiminus2=psiminus2*e2idelta+conj(psiminus2);
00399
00400 ddeta_psi=(psi2-psi0)/deleta;
00401 ddeta_psiminus=(psiminus2-psiminus0)/deleta;
00402
00403 I=(conj(psi)*psi+conj(psiminus)*psiminus)*a2*x/double(ell*ell);
00404 I-=conj(psiminus)*psi
00405 *((2.0*ell+1.0)*a2*double(ell)+2.0*eta*x*a2
00406 -eta*eta*double(ell))/(double(ell*ell)*root);
00407 I+=(conj(ddeta_psiminus)*psi-conj(ddeta_psi)*psiminus)
00408 *eta*root/double(ell);
00409
00410 }
00411 I=0.25*I;
00412 }
00413 else{
00414 x=q*epsilon/HBARC;
00415 if(ell==0){
00416 psi=sin(x+delta);
00417 psiplus=-cos(x+delta)+sin(x+delta)/x;
00418 psi=0.5*x*(Bessel::hstarn(0,x)+Misc::ceiphi(2.0*delta)*Bessel::hn(0,x));
00419 psiplus=0.5*x*(Bessel::hstarn(1,x)
00420 +Misc::ceiphi(2.0*delta)*Bessel::hn(1,x));
00421 psi*=Misc::ceiphi(-delta); psiplus*=Misc::ceiphi(-delta);
00422
00423
00424
00425
00426 I=(conj(psi)*psi+conj(psiplus)*psiplus)*x;
00427 I-=conj(psi)*psiplus;
00428 }
00429 else{
00430 psi=x*Bessel::hn(ell,x);
00431 psi=(Misc::ceiphi(2.0*delta)*psi+conj(psi))*0.5;
00432 psiminus=x*Bessel::hn(ell-1,x);
00433 psiminus=0.5*(Misc::ceiphi(2.0*delta)*psiminus+conj(psiminus));
00434
00435 I=(conj(psi)*psi+conj(psiminus)*psiminus)*x;
00436 I-=(2.0*double(ell)+1)*conj(psiminus)*psi;
00437
00438 }
00439 }
00440 return -real(I)/q;
00441 }
00442
00443 #endif