Particles and Nuclei Sub-Library: Version 0.910
frdm_mass.h
00001 /*
00002   -------------------------------------------------------------------
00003   
00004   Copyright (C) 2006-2012, Andrew W. Steiner
00005   
00006   This file is part of O2scl.
00007   
00008   O2scl is free software; you can redistribute it and/or modify
00009   it under the terms of the GNU General Public License as published by
00010   the Free Software Foundation; either version 3 of the License, or
00011   (at your option) any later version.
00012   
00013   O2scl is distributed in the hope that it will be useful,
00014   but WITHOUT ANY WARRANTY; without even the implied warranty of
00015   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00016   GNU General Public License for more details.
00017   
00018   You should have received a copy of the GNU General Public License
00019   along with O2scl. If not, see <http://www.gnu.org/licenses/>.
00020 
00021   -------------------------------------------------------------------
00022 */
00023 #ifndef FRDM_MASS_H
00024 #define FRDM_MASS_H
00025 
00026 #include <cmath>
00027 
00028 #include <o2scl/nucleus.h>
00029 #include <o2scl/nuclear_mass.h>
00030 #include <o2scl/constants.h>
00031 
00032 #ifndef DOXYGENP
00033 namespace o2scl {
00034 #endif
00035   
00036   /** \brief FRDM semi-empirical mass formula (macroscopic part only
00037       with no deformation)
00038       
00039       The macroscopic part of the finite-range droplet model
00040       from \ref Moller95 .
00041   
00042       Using the relations
00043       \f[
00044       \bar{\delta} = (n_n - n_p)/n
00045       \f]
00046       and
00047       \f[
00048       \bar{\epsilon} = - (n-n_0)/3/n_0
00049       \f]
00050       we get
00051       \f[
00052       n_n = \frac{1}{2} (1+\bar{\delta}) (1-3 \bar{\epsilon}) n_0
00053       \f]
00054       and
00055       \f[
00056       n_p = \frac{1}{2} (1-\bar{\delta}) (1-3 \bar{\epsilon}) n_0
00057       \f]
00058       Assuming that 
00059       \f[
00060       \frac{4 \pi}{3} R_n^3 n_n = N
00061       \f]
00062       and 
00063       \f[
00064       \frac{4 \pi}{3} R_p^3 n_p = Z
00065       \f]
00066       we get
00067       \f[
00068       R_n^3 = 3 N / \alpha_n
00069       \f]
00070       \f[
00071       R_p^3 = 3 Z / \alpha_p
00072       \f]
00073       where \f$ \alpha \f$'s are
00074       \f[
00075       \alpha_n = 2 \pi (1+ \bar{\delta})(1 - 3 \bar{\epsilon}) n_0
00076       \f]
00077       \f[
00078       \alpha_p = 2 \pi (1- \bar{\delta})(1 - 3 \bar{\epsilon}) n_0
00079       \f]
00080       Note that the above relations are somehow self-consistent
00081       because they imply
00082       \f[
00083       R^3 n = R_n^3 n_n + R_p^3 n_p
00084       \f]
00085 
00086       Since we're using (is there a better way?)
00087       \f[
00088       R = r_0 A^{1/3}
00089       \f]
00090       with \f$ r_0 = 1.16 \f$ fm, then 
00091       \f$ n_0 = 0.152946 \mathrm{fm}^{-3} \f$.
00092 
00093       \todo Fix pairing energy and double vs. int
00094       \todo Document drip_binding_energy(), etc.
00095       \todo Decide on number of fit parameters (10 or 12?) or
00096       let the user decide
00097       \todo Document the protected variables
00098       \todo Set the neutron and proton masses and hbarc to Moller et al.'s 
00099       values
00100 
00101       \future Add microscopic part.
00102   */
00103   class frdm_mass : public nuclear_mass_fit {
00104   public:
00105 
00106     /// Volume-energy constant in MeV (default 16.247)
00107     double a1;
00108     /// Symmetry-energy constant in MeV (default 32.73)
00109     double J;
00110     /// Nuclear compressibility constant in MeV (default 240)
00111     double K;
00112     /// Surface-energy constant in MeV (default 22.92)
00113     double a2;
00114     /// Effective surface-stiffness constant in MeV (default 29.21)
00115     double Q;
00116     /// Curvature-energy constant in MeV (default 0)
00117     double a3;
00118     /// Charge-asymmetry constant in MeV (default 0.436)
00119     double ca;
00120     /// Wigner constant in MeV (default 30)
00121     double W;
00122     /** \brief electronic-binding constant in MeV 
00123         (default \f$ 1.433 \times 10^{-5} \f$ ).
00124     */
00125     double ael;
00126     /// Proton root-mean-square radius in fm (default 0.80)
00127     double rp;
00128     /// Nuclear-radius constant in fm (default 1.16)
00129     double r0;
00130     /// Hydrogen atom mass excess, 7.289034 MeV
00131     double MH;
00132     /// Neutron mass excess, 8.071431 MeV
00133     double Mn;
00134     /// Electronic charge squared, 1.4399764 MeV fm
00135     double e2;
00136     /// Range of Yukawa-plus-exponential potential, 0.68 fm
00137     double a;
00138     /** \brief Range of Yukawa function used to generate nuclear 
00139         charge distribution, 0.70 fm
00140     */
00141     double aden;
00142     /// Average pairing-gap constant, 4.80 MeV
00143     double rmac;
00144     /// Neutron-proton interaction constant, 6.6 MeV
00145     double h;
00146     /// Density-symmetry constant, 0 MeV
00147     double L;
00148     /// Pre-exponential compressibility-term constant, 60 MeV
00149     double C;
00150     /// Exponential compressibility-term range constant, 0.831
00151     double gamma;
00152     /// Atomic mass unit, 931.5014 MeV
00153     double amu;
00154     
00155     /// Internal average neutron density
00156     double nn;
00157     
00158     /// Internal average proton density
00159     double np;
00160 
00161     /// Neutron radius
00162     double Rn;
00163 
00164     /// Proton radius
00165     double Rp;
00166 
00167     frdm_mass() {
00168       MH=7.289034;
00169       Mn=8.071431;
00170       e2=1.4399764;
00171       ael=1.433e-5;
00172       K=240.0;
00173       rp=0.8;
00174       r0=1.16;
00175       a=0.68;
00176       aden=0.7;
00177       rmac=4.8;
00178       h=6.6;
00179       W=30.0;
00180       L=0;
00181       a3=0;
00182       a1=16.247;
00183       a2=22.92;
00184       J=32.73;
00185       Q=29.21;
00186       a0=0.0;
00187       ca=0.436;
00188       C=60.0;
00189       gamma=0.831;
00190       amu=931.5014;
00191       //nfit=12;
00192       nfit=10;
00193     }
00194 
00195     /// Given \c Z and \c N, return the mass excess in MeV
00196     virtual double mass_excess_d(double Z, double N) {
00197       double ret;
00198       
00199       double A=Z+N;
00200       double cN=cbrt(N);
00201       double cZ=cbrt(Z);
00202       double cA=cbrt(A);
00203       double I=(N-Z)/A;
00204       
00205       // Ansatz for computing nuclear radius
00206       // This should imply a saturation density
00207       double R=r0*cA;
00208 
00209       double y=R/aden;
00210       double y0=r0*cA/aden;
00211       double x=R/a;
00212       double x0=r0*cA/a;
00213 
00214       double x02=x0*x0;
00215       B1=1.0-3.0/x02+(1.0+x0)*(2.0+3.0/x0+3.0/x02)*exp(-2.0*x0);
00216       B2=1.0-(1.0+2.0*x0-x02)*exp(-2.0*x0);
00217       double y02=y0*y0, y03=y02*y0, y04=y03*y0, y05=y04*y0;
00218       B3=1.0-5.0/y02*(1.0-15.0/8.0/y0+21.0/8.0/y03-0.75*
00219                       (1.0+4.5/y0+7.0/y02+3.5/y03)*exp(-2.0*y0));
00220       B4=1.0+5.0*(-3.0/y02+7.5/y03-63.0/4.0/pow(y0,5.0)
00221                   +0.75*(2.0/y0+12.0/y02+32.0/y03+
00222                          42.0/y04+21.0/y05)*exp(-2.0*y0));
00223       
00224       // Assume spherical shape
00225       Bs=1.0;
00226       Bv=1.0;
00227       Bw=1.0;
00228       Bk=1.0;
00229       Br=1.0;
00230       
00231       /* Pairing parameters:
00232 
00233          Deltan goes like N^{-1/3}
00234          Deltap goes like Z^{-1/3}
00235          deltanp goes like A^{-2/3}
00236       */
00237       Deltan=rmac*Bs/cN;
00238       Deltap=rmac*Bs/cZ;
00239       deltanp=h/Bs/cA/cA;
00240 
00241       // Coulomb coefficient
00242       c1=0.6*e2/r0;
00243       // Volume redistribution coefficient
00244       c2=1.0/336.0*(1.0/J+18.0/K)*c1*c1;
00245       // Coulomb exchange coefficient
00246       c4=1.25*pow(3.0/2.0/o2scl_const::pi,2.0/3.0)*c1;
00247       // Surface redistribution coefficient
00248       c5=c1*c1/64.0/Q;
00249       // Proton form-factor coefficient
00250       f0=-0.125*145.0/48.0*rp*rp*e2/r0/r0/r0;
00251 
00252       // Average bulk nuclear asymmetry, which
00253       // goes like \f$ (I,Z A^{-2/3} B1^{-1}) / (1,B1^{-1}) \f$.
00254       deltabar=(I+3.0/16.0*c1/Q*Z/cA/cA*Bv*Bs/B1)/
00255         (1.0+2.25*J/Q/cA*Bs*Bs/B1);
00256 
00257       // Average relative deviation of bulk density
00258       epsbar=(C*exp(-gamma*cA)-2.0*a2*B2/cA+L*deltabar*deltabar+
00259               c1*Z*Z/A/cA*B4)/K;
00260       
00261       // Pairing contribution
00262       double tm=0.0, tm2;
00263       if (((int)Z)%2==1 && ((int)N)%2==1) {
00264         if (Z==N) tm=1.0/A;
00265         tm2=Deltap+Deltan-deltanp;
00266       } else if (((int)Z)%2==1 && ((int)N)%2==0) {
00267         tm2=Deltap;
00268       } else if (((int)N)%2==1) {
00269         tm2=Deltan;
00270       } else {
00271         tm2=0.0;
00272       }
00273       
00274       ret=MH*Z+Mn*N+(-a1+J*deltabar*deltabar-0.5*K*epsbar*epsbar)*A
00275         +(a2*B1+2.25*J*J/Q*deltabar*deltabar*Bs*Bs/B1)*cA*cA
00276         +a3*cA*Bk+c1*Z*Z/cA*B3-c2*Z*Z*cA*Br-c4*Z*cZ/cA
00277         -c5*Z*Z*Bw*Bs/B1+f0*Z*Z/A-ca*(N-Z)+W*(fabs(I)+tm)+tm2
00278         -ael*pow(Z,2.39);
00279       
00280       return ret;
00281     }
00282 
00283     virtual int fit_fun(size_t nv, const ovector_base &x) {
00284       K=x[0]*200.0;
00285       r0=x[1];
00286       W=x[2];
00287       a1=x[3];
00288       a2=x[4];
00289       J=x[5];
00290       Q=x[6];
00291       ca=x[7];
00292       C=x[8]*60.0;
00293       gamma=x[9];
00294       //rmac=x[10];
00295       //h=x[11];
00296       return 0;
00297 
00298     }
00299 
00300     virtual int guess_fun(size_t nv, ovector_base &x) {
00301       x[0]=K/200.0;
00302       x[1]=r0;
00303       x[2]=W;
00304       x[3]=a1;
00305       x[4]=a2;
00306       x[5]=J;
00307       x[6]=Q;
00308       x[7]=ca;
00309       x[8]=C/60.0;
00310       x[9]=gamma;
00311       //x[10]=rmac;
00312       //x[11]=h;
00313       return 0;
00314 
00315     }
00316 
00317     /** \brief Return the binding energy in MeV
00318      */
00319     virtual double drip_binding_energy_d(double Z, double N,
00320                                          double npout, double nnout,
00321                                          double chi) {
00322       return (drip_mass_excess_d(Z,N,npout,nnout,chi)+
00323               ((Z+N)*o2scl_fm::mass_amu-Z*o2scl_fm::mass_electron-
00324                N*o2scl_fm::mass_neutron-Z*o2scl_fm::mass_proton)*
00325               o2scl_const::hc_mev_fm);
00326     }
00327 
00328     /** \brief Given \c Z and \c N, return the mass excess in MeV
00329      */
00330     virtual double drip_mass_excess_d(double Z, double N,
00331                                       double np_out, double nn_out,
00332                                       double chi) {
00333       double ret;
00334       
00335       double dN=N;
00336       double dZ=Z;
00337       double dA=Z+N;
00338       double cN=cbrt(dN);
00339       double cZ=cbrt(dZ);
00340       double cA=cbrt(dA);
00341       double I=(dN-dZ)/dA;
00342       
00343       // Ansatz for computing nuclear radius
00344       // This should imply a saturation density
00345       double R=r0*cA;
00346 
00347       double y=R/aden;
00348       double y0=r0*cA/aden;
00349       double x=R/a;
00350       double x0=r0*cA/a;
00351 
00352       double x02=x0*x0;
00353       B1=1.0-3.0/x02+(1.0+x0)*(2.0+3.0/x0+3.0/x02)*exp(-2.0*x0);
00354       B2=1.0-(1.0+2.0*x0-x02)*exp(-2.0*x0);
00355       double y02=y0*y0, y03=y02*y0, y04=y03*y0, y05=y04*y0;
00356       B3=1.0-5.0/y02*(1.0-15.0/8.0/y0+21.0/8.0/y03-0.75*
00357                       (1.0+4.5/y0+7.0/y02+3.5/y03)*exp(-2.0*y0));
00358       B4=1.0+5.0*(-3.0/y02+7.5/y03-63.0/4.0/pow(y0,5.0)
00359                   +0.75*(2.0/y0+12.0/y02+32.0/y03+
00360                          42.0/y04+21.0/y05)*exp(-2.0*y0));
00361       
00362       // Assume spherical shape
00363       Bs=1.0;
00364       Bv=1.0;
00365       Bw=1.0;
00366       Bk=1.0;
00367       Br=1.0;
00368       
00369       /* Pairing parameters:
00370 
00371          Deltan goes like N^{-1/3}
00372          Deltap goes like Z^{-1/3}
00373          deltanp goes like A^{-2/3}
00374       */
00375       Deltan=rmac*Bs/cN;
00376       Deltap=rmac*Bs/cZ;
00377       deltanp=h/Bs/cA/cA;
00378 
00379       // Coulomb coefficient
00380       c1=0.6*e2/r0;
00381       // Volume redistribution coefficient
00382       c2=1.0/336.0*(1.0/J+18.0/K)*c1*c1;
00383       // Coulomb exchange coefficient
00384       c4=1.25*pow(3.0/2.0/o2scl_const::pi,2.0/3.0)*c1;
00385       // Surface redistribution coefficient
00386       c5=c1*c1/64.0/Q;
00387       // Proton form-factor coefficient
00388       f0=-0.125*145.0/48.0*rp*rp*e2/r0/r0/r0;
00389 
00390       // Average bulk nuclear asymmetry, which
00391       // goes like \f$ (I,Z A^{-2/3} B1^{-1}) / (1,B1^{-1}) \f$.
00392       deltabar=(I+3.0/16.0*c1/Q*dZ/cA/cA*Bv*Bs/B1)/
00393         (1.0+2.25*J/Q/cA*Bs*Bs/B1);
00394 
00395       // Average relative deviation of bulk density
00396       epsbar=(C*exp(-gamma*cA)-2.0*a2*B2/cA+L*deltabar*deltabar+
00397               c1*dZ*dZ/dA/cA*B4)/K;
00398       
00399       // Saturation density
00400       double n0=0.152946;
00401 
00402       // Determine densities and radii
00403       nn=0.5*(1.0+deltabar)*(1.0-3.0*epsbar)*n0;
00404       np=0.5*(1.0-deltabar)*(1.0-3.0*epsbar)*n0;
00405       
00406       Rn=cbrt(3*N/2.0/o2scl_const::pi/(1.0+deltabar)/(1.0-3.0*epsbar)/n0);
00407       Rp=cbrt(3*Z/2.0/o2scl_const::pi/(1.0-deltabar)/(1.0-3.0*epsbar)/n0);
00408       
00409       double tm=0.0, tm2=0.0;
00410       if (false) {
00411         // Pairing contribution
00412         if (((int)Z)%2==1 && ((int)N)%2==1) {
00413           if (Z==N) tm=1.0/dA;
00414           tm2=Deltap+Deltan-deltanp;
00415         } else if (((int)Z)%2==1 && ((int)N)%2==0) {
00416           tm2=Deltap;
00417         } else if (((int)N)%2==1) {
00418           tm2=Deltan;
00419         } else {
00420           tm2=0.0;
00421         }
00422       }
00423       
00424       //double xp=np/(nn+np);
00425       double x3fact=1.0;//+(xp*xp*xp-1.0)/(1.0+exp((xp-0.25)/0.01));
00426       
00427       double surf_part=(a2*B1+2.25*J*J/Q*deltabar*deltabar*Bs*Bs/B1)*cA*cA+
00428         -c5*dZ*dZ*Bw*Bs/B1;
00429       
00430       ret=MH*dZ+Mn*dN+(-a1+J*deltabar*deltabar-0.5*K*epsbar*epsbar)*dA
00431         +a3*cA*Bk+c1*dZ*dZ/cA*B3-c2*dZ*dZ*cA*Br-c4*dZ*cZ/cA
00432         +f0*dZ*dZ/dA-ca*(dN-dZ)+W*(fabs(I)+tm)+tm2
00433         -ael*pow(Z,2.39)+surf_part*x3fact;
00434       
00435       return ret;
00436     }
00437 
00438   protected:
00439 
00440     /// Proton pairing coefficient
00441     double Deltap;
00442     /// Neutron pairing coefficient
00443     double Deltan;
00444     /// Isovector pairing coefficient
00445     double deltanp;
00446 
00447     /// Average bulk nuclear asymmetry
00448     double deltabar;
00449     /// Average relative deviation of bulk density
00450     double epsbar;
00451 
00452     /// Desc
00453     double Bs;
00454     /// Desc
00455     double Bk;
00456     /// Desc
00457     double Br;
00458     /// Desc
00459     double Bw;
00460     /// Desc
00461     double Bv;
00462 
00463     /// Coulomb energy coefficient
00464     double c1;
00465     /// Volume redistribution energy coefficient
00466     double c2;
00467     /// Coulomb exchange correction coefficient
00468     double c4;
00469     /// Surface redistribution energy coefficient
00470     double c5;
00471 
00472     /** \brief Coefficient for the proton form-factor correction to the 
00473         Coulomb energy
00474     */
00475     double f0;
00476 
00477     /// Desc
00478     double a0;
00479 
00480     /// Desc
00481     double B1;
00482     /// Desc
00483     double B2;
00484     /// Desc
00485     double B3;
00486     /// Desc
00487     double B4;
00488 
00489   };
00490   
00491 #ifndef DOXYGENP
00492 }
00493 #endif
00494 
00495 #endif
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Friends

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).

Get Object-oriented Scientific Computing
Lib at SourceForge.net. Fast, secure and Free Open Source software
downloads.