All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
nucleus_rmf.h
Go to the documentation of this file.
1 /*
2  -------------------------------------------------------------------
3 
4  Copyright (C) 2006-2014, Andrew W. Steiner
5 
6  This file is part of O2scl.
7 
8  O2scl is free software; you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version.
12 
13  O2scl is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with O2scl. If not, see <http://www.gnu.org/licenses/>.
20 
21  -------------------------------------------------------------------
22 */
23 /** \file nucleus_rmf.h
24  \brief File defining \ref o2scl::nucleus_rmf
25 */
26 #ifndef RMF_NUCLEUS_H
27 #define RMF_NUCLEUS_H
28 
29 #include <iostream>
30 #include <string>
31 #include <vector>
32 #include <o2scl/interp.h>
33 #include <o2scl/constants.h>
34 #include <o2scl/part.h>
35 #include <o2scl/eos_had_rmf.h>
36 #include <o2scl/table_units.h>
37 #include <o2scl/ode_rkck_gsl.h>
38 #include <o2scl/ode_funct.h>
39 #include <o2scl/shared_ptr.h>
40 
41 #ifndef DOXYGEN_NO_O2NS
42 namespace o2scl {
43 #endif
44 
45  /** \brief Spherical closed-shell nuclei with a relativistic
46  mean-field model in the Hartree approximation
47 
48  This code is very experimental.
49 
50  This class is based on a code developed by C.J. Horowitz and
51  B.D. Serot, and used in \ref Horowitz81 which was then adapted
52  by P.J. Ellis and used in \ref Heide94 and \ref Prakash94. Ellis
53  and A.W. Steiner adapted it for the parameterization in in \ref
54  eos_had_rmf for \ref Steiner05b, and then converted to C++ by
55  Steiner afterwards.
56 
57  The standard usage is something like:
58  \code
59  nucleus_rmf rn;
60  o2scl_hdf::rmf_load(rn.rmf,"NL4");
61  rn.run_nucleus(82,208,0,0);
62  cout << rn.rnrp << endl;
63  \endcode
64  which computes the structure of \f$ ^{208}\mathrm{Pb} \f$ and
65  outputs the neutron skin thickness using the model \c 'NL4'.
66 
67  Potential exceptions are
68  - Failed to converge
69  - Failed to solve meson field equations
70  - Energy not finite (usually a problem in the equation of
71  state)
72  - Energy not finite in final calculation
73  - Function \ref iterate() called before \ref init_run()
74  - Not a closed-shell nucleus
75 
76  The initial level pattern is
77  \verbatim
78  1 S 1/2
79  // 2 nucleons
80  1 P 3/2
81  1 P 1/2
82  // 8 nucleus
83  1 D 5/2
84  1 D 3/2
85  2 S 1/2
86  // 20 nucleons
87  1 F 7/2
88  // 28 nucleons
89  1 F 5/2
90  2 P 3/2
91  2 P 1/2
92  // 40 nucleons
93  1 G 9/2
94  // 50 nucleus
95  1 G 7/2
96  2 D 5/2
97  1 H 11/2
98  2 D 3/2
99  3 S 1/2
100  // 82 nucleons
101  1 H 9/2
102  2 F 7/2
103  1 I 13/2
104  2 F 5/2
105  3 P 3/2
106  3 P 1/2
107  // 126 nucleons
108  2 G 9/2
109  1 I 11/2
110  1 J 15/2
111  3 D 5/2
112  4 S 1/2
113  2 G 7/2
114  3 D 3/2
115  // 184 nucleons
116  \endverbatim
117 
118  Below, \f$ \alpha \f$ is a generic index for the isospin, the
119  radial quantum number \f$ n \f$ and the angular quantum numbers
120  \f$ \kappa \f$ and \f$ m \f$. The meson fields are \f$
121  \sigma(r), \omega(r) \f$ and \f$ \rho(r) \f$. The baryon density
122  is \f$ n(r) \f$, the neutron and proton densities are \f$ n_n(r)
123  \f$ and \f$ n_p(r) \f$, and the baryon scalar density is \f$
124  n_s(r) \f$.
125  The nucleon field equations are
126  \f{eqnarray*}
127  F^{\prime}_{\alpha}(r)- \frac{\kappa}{r} F_{\alpha}(r)
128  + \left[ \varepsilon_{\alpha} - g_{\omega} \omega(r)
129  - t_{\alpha} g_{\rho} \rho(r) - \left( t_{\alpha}+\frac{1}{2} \right)
130  e A(r) - M + g_{\sigma} \sigma(r) \right] G_{\alpha}(r) &=& 0 \\
131  G^{\prime}_{\alpha}(r)+ \frac{\kappa}{r} G_{\alpha}(r)
132  - \left[ \varepsilon_{\alpha} - g_{\omega} \omega(r)
133  - t_{\alpha} g_{\rho} \rho(r) - \left( t_{\alpha}+\frac{1}{2}
134  \right) e A(r) + M - g_{\sigma} \sigma(r) \right] F_{\alpha}(r) &=& 0
135  \f}
136  where the isospin number, \f$ t_{\alpha} \f$ is \f$ 1/2 \f$ for
137  protons and \f$ -1/2 \f$ for neutrons.
138  The meson field equations are
139  \f{eqnarray*}
140  \sigma^{\prime \prime}(r) + \frac{2}{r} \sigma^{\prime}(r)
141  - m_{\sigma}^2 \sigma &=& - g_{\sigma} n_s(r) + b M g_{\sigma}^3
142  \sigma^2 + c g_{\sigma}^4 \sigma^3 - g_{\rho}^2 \rho^2
143  \frac{\partial f}{\partial \sigma} \\
144  \omega^{\prime \prime}(r) + \frac{2}{r} \omega^{\prime}(r)
145  - m_{\omega}^2 \omega &=& - g_{\omega} n(r) +
146  \frac{\zeta}{6} g_{\omega}^4 \omega^3 + g_{\rho}^2 \rho^2
147  \frac{\partial f}{\partial \omega} \\
148  \rho^{\prime \prime}(r) + \frac{2}{r} \rho^{\prime}(r)
149  - m_{\rho}^2 \rho &=& - \frac{g_{\rho}}{2}
150  \left[n_n(r)-n_p(r)\right] + 2 g_{\rho}^2 \rho f + \frac{\xi}{6}
151  g_{\rho}^4 \rho^3
152  \f}
153  and the Coulomb field equation is
154  \f[
155  A^{\prime \prime}(r) + \frac{2}{r} A^{\prime}(r) = - e n_p(r)
156  \f]
157  The meson field equations plus a pair of Dirac-like nucleon
158  field equations for each index \f$ \alpha \f$ must be solved to
159  compute the structure of a given nucleus.
160 
161  The densities (scalar, baryon, isospin, and charge) are
162  \f{eqnarray*}
163  n_s(r) &=& \sum_{\alpha} \left\{ \int d^3 r \left[ G(r)^2-F(r)^2
164  \right] \right\} \\
165  n(r) &=& \sum_{\alpha} \left\{ \int d^3 r \left[ G(r)^2+F(r)^2
166  \right] \right\} \\
167  n_i(r) &=& \sum_{\alpha} \left\{ t_{\alpha} \int d^3 r
168  \left[ G(r)^2-F(r)^2 \right] \right\} \\
169  n_c(r) &=& \sum_{\alpha} \left\{ \left[t_{\alpha}+\frac{1}{2}\right]
170  \int d^3 r \left[ G(r)^2-F(r)^2 \right] \right\}
171  \f}
172 
173  The total energy is
174  \f[
175  E = \sum_{\alpha} \varepsilon_{\alpha} (2 j_{\alpha}+1)
176  - \frac{1}{2} \int d^{3} r
177  \left[- g_{\sigma} \sigma(r) \rho_s(r) + g_{\omega} \omega(r)
178  \rho(r) + \frac{1}{2} g_{\rho} \rho(r) + e A(r) n_p(r) \right]
179  \f]
180 
181  The charge density is the proton density folded with the
182  charge density of the proton, i.e.
183  \f[
184  \rho_{\mathrm{ch}}(r) = \int d^{3} r^{\prime}
185  \rho_{\mathrm{prot}}(r-r^{\prime}) \rho_p(r)
186  \f]
187  where the proton charge density is assumed to be of the form
188  \f[
189  \rho_{\mathrm{prot}}(r) = \frac{\mu^3}{8 \pi} \exp \left(
190  - \mu |r|\right)
191  \f]
192  and the parameter \f$ \mu = (0.71)^{1/2}~\mathrm{GeV} \f$ (see
193  Eq. 20b in \ref Horowitz81). The default value of \ref a_proton
194  is the value of \f$ \mu \f$ converted into \f$ \mathrm{fm}^{-1}
195  \f$.
196 
197  Generally, the first array index associated with a function
198  is not the value at \f$ r=0 \f$, but at \f$ r=\Delta r \f$
199  (stored in \ref step_size) and so the \f$ r=0 \f$ part of the
200  algorithm is handled separately.
201 
202  \todo Better documentation
203  \todo Convert energies() to use EOS and possibly
204  replace sigma_rhs() and related functions by the associated
205  field equation method of eos_had_rmf.
206 
207  \todo Document hw=3.923+23.265/cbrt(atot);
208 
209  \comment
210  the hw=blah blah term for the CM correction is discussed
211  a bit in Negele, PRC 1 (1970) 1260, esp. see Eq. 2.30 but
212  the numerical coefficients are different here.
213  \endcomment
214 
215  \todo I believe currently the center of mass correction
216  for the binding energy per nucleon is not done and has
217  to be added in after the fact
218 
219  \future Sort energy levels at the end by energy
220  \future Improve the numerical methods
221  \future Make the neutron and proton orbitals more configurable
222  \future Generalize to \f$ m_n \neq m_p \f$ .
223  \future Allow more freedom in the integrations
224  \future Consider converting everything to inverse fermis.
225  \future Convert to zero-indexed arrays (mostly done)
226  \future Warn when the level ordering is wrong, and unoccupied levels
227  are lower energy than occupied levels
228  \future Connect with \ref o2scl::nucmass ?
229  */
230  class nucleus_rmf {
231 
234 
235 #ifndef DOXYGEN_INTERNAL
236 
237  protected:
238 
239  /// A convenient struct for the solution of the Dirac equations
240  typedef struct {
241  /// Eigenvalue
242  double eigen;
243  /// Quantum number \f$ \kappa \f$ .
244  double kappa;
245  /// The meson fields
247  /// Desc
249  } odparms;
250 
251  /// The total number of shells stored internally
252  static const int n_internal_levels=29;
253 
254 #endif
255 
256  public:
257 
258  nucleus_rmf();
259 
260  ~nucleus_rmf();
261 
262  /** \name Basic operation
263  */
264  //@{
265  /// A shell of nucleons for \ref nucleus_rmf
266  typedef struct {
267  /// Degeneracy \f$ 2 j+1 \f$ .
268  int twojp1;
269  /// \f$ \kappa \f$
270  int kappa;
271  /// Energy eigenvalue
272  double energy;
273  /// Isospin ( \f$ +1/2 \f$ or \f$ -1/2 \f$ .
274  double isospin;
275  /// Angular momentum-spin state \f$ ^{2s+1} \ell_{j} \f$
276  std::string state;
277  /// Matching radius (in fm)
278  double match_point;
279  /// Desc.
280  double eigen;
281  /// Desc.
282  double eigenc;
283  /// Number of nodes in the wave function
284  int nodes;
285  } shell;
286 
287  /** \brief Computes the structure of a nucleus with the specified
288  number of levels
289 
290  Note that \ref rmf must be set before run_nucleus() is called.
291 
292  This calls init_run(), and then iterate() until \c iconverged is
293  1, and then post_converge().
294  */
295  int run_nucleus(int nucleus_Z, int nucleus_N, int unocc_Z, int unocc_N);
296 
297  /// Set output level
298  void set_verbose(int v) { verbose=v; };
299  //@}
300 
301  /** \name Lower-level interface
302  */
303  //@{
304  /** \brief Initialize a run
305 
306  Note that \ref rmf must be set before run_nucleus() is called.
307  */
308  void init_run(int nucleus_Z, int nucleus_N, int unocc_Z, int unocc_N);
309 
310  /// Perform an iteration
311  void iterate(int nucleus_Z, int nucleus_N, int unocc_Z, int unocc_N,
312  int &iconverged);
313 
314  /// After convergence, make CM corrections, etc.
315  int post_converge(int nucleus_Z, int nucleus_N, int unocc_Z, int unocc_N);
316 
317  //@}
318 
319  /// \name Results
320  //@{
321  /** \brief Get the radial profiles
322 
323  The profiles are calculated each iteration by iterate().
324  */
326 
327  /** \brief The final charge densities
328  */
330 
331  /// The number of levels
332  int nlevels;
333 
334  /** \brief The levels (protons first, then neutrons)
335 
336  An array of size \ref nlevels
337  */
338  std::vector<shell> levels;
339 
340  /// The number of unoccupied levels (equal to \c unocc_Z + \c unocc_N)
342 
343  /** \brief The unoccupied levels (protons first, then neutrons)
344 
345  An array of size \ref nuolevels
346  */
347  std::vector<shell> unocc_levels;
348 
349  /** \brief Surface tension (in \f$ \mathrm{fm}^{-3} \f$ )
350 
351  Computed in post_converge() or automatically in run_nucleus()
352  */
353  double stens;
354 
355  /** \brief Skin thickness (in fm)
356 
357  Computed every iteration in iterate()
358  or automatically in run_nucleus()
359  */
360  double rnrp;
361 
362  /** \brief Neutron RMS radius (in fm)
363 
364  Computed every iteration in iterate() or automatically in
365  run_nucleus()
366  */
367  double rnrms;
368 
369  /** \brief Proton RMS radius (in fm)
370 
371  Computed every iteration in iterate() or automatically in
372  run_nucleus()
373  */
374  double rprms;
375 
376  /** \brief Total energy (in MeV)
377 
378  Computed every iteration in iterate() or automatically in
379  run_nucleus()
380  */
381  double etot;
382 
383  /** \brief Charge radius (in fm)
384 
385  Computed in post_converge() or automatically in run_nucleus()
386  */
387  double r_charge;
388 
389  /** \brief Charge radius corrected by the center of mass (in fm)
390 
391  Computed in post_converge() or automatically in run_nucleus()
392  */
393  double r_charge_cm;
394  //@}
395 
396  /** \name Equation of state
397  */
398  //@{
399  /** \brief The default equation of state (default NL3)
400 
401  This is set in the constructor to be the default
402  model, NL3, using the function \ref load_nl3().
403  */
405 
406  /** \brief Set the base EOS to be used
407 
408  The equation of state must be set before run_nucleus() or
409  init_fun() are called, including the value of eos_had_rmf::mnuc.
410  */
412  rmf=&r;
413  return 0;
414  }
415 
416  /** \brief \ref thermo object for the EOS
417 
418  This is just used as temporary storage.
419  */
421 
422  /** \brief The neutron
423 
424  The mass of the neutron is ignored and set by init_run()
425  to be eos_had_rmf::mnuc from \ref rmf.
426  */
428 
429  /** \brief The proton
430 
431  The mass of the proton is ignored and set by init_run()
432  to be eos_had_rmf::mnuc from \ref rmf.
433  */
435  //@}
436 
437  /** \name Numeric configuration
438  */
439  //@{
440  /** \brief If true, call the error handler if the routine does not
441  converge or reach the desired tolerance (default true)
442  */
444 
445  /// Set the stepper for the Dirac differential equation
447  ode_funct11> &step) {
448  ostep=&step;
449  return;
450  }
451 
452  /// Maximum number of total iterations (default 70)
453  int itmax;
454 
455  /** \brief Maximum number of iterations for solving the meson
456  field equations (default 10000)
457  */
459 
460  /** \brief Maximum number of iterations for solving the Dirac
461  equations (default 100)
462  */
464 
465  /// Tolerance for Dirac equations (default \f$ 5 \times 10^{-3} \f$ ).
466  double dirac_tol;
467 
468  /** \brief Second tolerance for Dirac equations
469  (default \f$ 5 \times 10^{-4} \f$ ).
470  */
471  double dirac_tol2;
472 
473  /// Tolerance for meson field equations (default \f$ 10^{-6} \f$ ).
474  double meson_tol;
475  //@}
476 
477  /** \brief Initial guess structure
478 
479  The initial guess for the meson field profiles is
480  a set of Fermi-Dirac functions, i.e.
481  \f[
482  \sigma(r)=\mathrm{sigma0}/
483  [1+\exp((r-\mathrm{fermi\_radius})/\mathrm{fermi\_width})]
484  \f]
485  */
486  typedef struct {
487  /// Scalar field at r=0
488  double sigma0;
489  /// Vector field at r=0
490  double omega0;
491  /// Isubvector field at r=0
492  double rho0;
493  /// Coulomb field at r=0
494  double A0;
495  /// The radius for which the fields are half their central value
496  double fermi_radius;
497  /// The "width" of the Fermi-Dirac function
498  double fermi_width;
499  } initial_guess;
500 
501  /** \brief Parameters for initial guess
502 
503  Default is {310,240,-6,25.9,6.85,0.6}
504  */
506 
507  /** \brief If true, use the generic ODE solver instead of the
508  internal 4th order Runge-Kutta
509  */
511 
512 #ifndef DOXYGEN_INTERNAL
513 
514  protected:
515 
516  /** \brief The parameter for the charge density of the proton
517  (default is about 4.27073)
518  */
519  double a_proton;
520 
521  /// Load the default model NL3 into the given \ref eos_had_rmf object
522  int load_nl3(eos_had_rmf &r);
523 
524  /// The base EOS
526 
527  /** \brief The radial profiles
528  */
530 
531  /** \brief The final charge densities
532  */
534 
535  /** \brief A pointer to the current vector of levels
536  (either \ref levels or \ref unocc_levels)
537  */
538  std::vector<shell> *levp;
539 
540  /** \brief Control output (default 1)
541  */
542  int verbose;
543 
544  /// The starting neutron levels
546 
547  /// The starting proton levels
549 
550  /// The grid size
551  static const int grid_size=300;
552 
553  /// The grid step size (default 0.04)
554  double step_size;
555 
556  /// The nucleon mass (automatically set in init_fun())
557  double mnuc;
558 
559  /** \name The meson and photon fields and field equations (protected)
560  */
561  //@{
562  /// Values of the fields from the last iteration
564 
565  /// The values of the fields
567 
568  /// The Green's functions inside
570 
571  /// The Green's functions outside
573 
574  /// Scalar density RHS
575  double sigma_rhs(double sig, double ome, double rho);
576 
577  /// Vector density RHS
578  double omega_rhs(double sig, double ome, double rho);
579 
580  /// Isubvector density RHS
581  double rho_rhs(double sig, double ome, double rho);
582 
583  /** \brief Calculate meson and photon
584  Green's functions \ref gin and \ref gout
585  */
586  void meson_init();
587 
588  /** \brief Calculate meson and photon fields
589 
590  The meson and photon field equations are of the
591  Sturm-Liouville form, e.g.
592  \f[
593  \left[r^2 \sigma^{\prime}(r) \right]^{\prime} -
594  r^2 m_{\sigma}^2 \sigma(r) = f(r)
595  \f]
596  where \f$ \sigma(0) = \sigma_0 \f$ and \f$ \sigma(+\infty) = 0
597  \f$. A solution of the homogeneous equation with \f$ f(r) =0
598  \f$ is \f$ \sigma(r) = \sigma_0 \sinh( m_{\sigma} r)/
599  (m_{\sigma} r) \f$. The associated Green's function is
600  \f[
601  D(r,r^{\prime},m_{\sigma}) = \frac{-1}{m_{\sigma} r r^{\prime}}
602  \sinh (m_{\sigma} r_{<}) \exp (-m_{\sigma} r_{>})
603  \f]
604  where \f$ r_{<} \f$ is the smaller of \f$ r \f$ and
605  \f$ r^{\prime} \f$ and \f$ r_{>} \f$ is the larger.
606  One can then write the solution of the meson field
607  equation given the density
608  \f[
609  \sigma(r) = \int_0^{\infty} r^{\prime 2}~d r^{\prime}
610  \left[ - g_{\sigma} n_{s}(r) \right]
611  D\left(r,r^{\prime},m_{\sigma}\right)
612  \f]
613 
614  Since radii are positive, \f$ \sinh (m r) \approx
615  e^{m r}/2 \f$ and
616  \f[
617  D(r,r^{\prime},m_{\sigma}) \approx \left[
618  \frac{-1}{2 m_{\sigma} r_{<}}
619  \exp (m_{\sigma} r_{<})
620  \right] \left[
621  \frac{1}{r_{>}}
622  \exp (-m_{\sigma} r_{>})
623  \right]
624  \f]
625  The two terms in the Green's function are separated into
626  \f[
627  \mathrm{gin}(r) = \frac{e^{m_{\sigma} r}}{2 m_{\sigma} r}
628  \f]
629  and
630  \f[
631  \mathrm{gout}(r) = \frac{e^{-m_{\sigma} r}}{r} \, .
632  \f]
633  These functions are computed in \ref meson_init() . Then the
634  field is given by
635  \f[
636  \sigma(r)= \left(\frac{e^{-m_{\sigma} r}}{r}\right)
637  \int_0^{r} r^{\prime 2} g_{\sigma} n_{s}
638  \left(\frac{e^{m_{\sigma} r^{\prime}}}{2 m_{\sigma}
639  r^{\prime}} \right)~d r^{\prime} +
640  \left(\frac{e^{m_{\sigma} r}}{2 m_{\sigma} r} \right)
641  \int_r^{\infty} r^{\prime 2} g_{\sigma} n_{s}
642  \left(\frac{e^{-m_{\sigma} r^{\prime}}}
643  {r^{\prime}}\right)~d r^{\prime}
644  \f]
645  where the first integral is stored in <tt>xi2</tt> and
646  the second is in <tt>xi1</tt> in the function \ref meson_iter() .
647  The part of <tt>xi2</tt> at \f$ r^{\prime}=0 \f$ is
648  stored in <tt>xi20</tt>.
649  */
650  void meson_iter(double ic);
651 
652  /** \brief Solve for the meson profiles
653  */
654  void meson_solve();
655 
656  /** \brief The grid index corresponding to the nuclear surface
657  (computed by init_run())
658  */
660  //@}
661 
662  /** \name Density information (protected)
663  */
664  //@{
665  /// The densities times radius squared
667 
668  /// The proton scalar density times radius squared
670 
671  /// The scalar field RHS
673 
674  /// The vector field RHS
676 
677  /// The isubvector field RHS
679 
680  /// Charge density
682 
683  /// Charge density
685 
686  /// Baryon density
688  //@}
689 
690  /// Energy integrand
692 
693  /// Initialize the meson and photon fields, the densities, etc.
694  void init_meson_density();
695 
696  /// Calculate the energy profile
697  void energies(double xpro, double xnu, double e);
698 
699  /// Compute the center of mass correction
700  void center_mass_corr(double atot);
701 
702  /** \name Calculating the form factor, etc. (protected)
703  */
704  //@{
705  /// Fold in proton form factor
706  void pfold(double x, double &xrhof);
707 
708  /// Function representing proton form factor
709  double xpform(double x, double xp, double a);
710 
711  /// Perform integrations for form factor
712  void gauss(double xmin, double xmax, double x, double &xi);
713 
714  /// Desc.
715  double xrhop(double x1);
716 
717  /// Interpolation object
719  //@}
720 
721  /** \name Solving the Dirac equations (protected)
722  */
723  //@{
724  /** \brief Solve the Dirac equations
725 
726  Solves the Dirac equation in from 12 fm to the match point and
727  then out from .04 fm and adjusts eigenvalue with
728  \f[
729  \Delta \varepsilon = -g(r=\mathrm{match\_point})
730  \times (f^{+}-f^{-})
731  \f]
732  */
733  void dirac(int ilevel);
734 
735  /// Take a step in the Dirac equations
736  void dirac_step(double &x, double h, double eigen,
737  double kappa, ubvector &varr);
738 
739  /// The form of the Dirac equations for the ODE stepper
740  int odefun(double x, size_t nv, const ubvector &y,
741  ubvector &dydx, odparms &op);
742 
743  /// Compute the fields for the Dirac equations
744  void field(double x, double &s, double &v, ubvector &varr);
745 
746  /// The default stepper
749 
750  /// The ODE stepper
753 
754  /** \brief Integrate the Dirac equations using a simple
755  inline 4th order Runge-Kutta
756  */
757  double dirac_rk4(double x, double g1, double f1, double &funt,
758  double eigen, double kappa, ubvector &varr);
759  //@}
760 
761  /// True if init() has been called
763 
764  /// ODE functions
766 
767  /// ODE derivatives
769 
770  /// ODE errors
772 
773  /// \name Gauss-Legendre integration points and weights
774  //@{
775  double x12[6], w12[6];
776  double x100[50], w100[50];
777  //@}
778 
779 #endif
780 
781  };
782 
783 #ifndef DOXYGEN_NO_O2NS
784 }
785 #endif
786 
787 #endif
void pfold(double x, double &xrhof)
Fold in proton form factor.
void dirac(int ilevel)
Solve the Dirac equations.
double dirac_tol
Tolerance for Dirac equations (default ).
Definition: nucleus_rmf.h:466
eos_had_rmf * rmf
The base EOS.
Definition: nucleus_rmf.h:525
double stens
Surface tension (in )
Definition: nucleus_rmf.h:353
int twojp1
Degeneracy .
Definition: nucleus_rmf.h:268
bool generic_ode
If true, use the generic ODE solver instead of the internal 4th order Runge-Kutta.
Definition: nucleus_rmf.h:510
A convenient struct for the solution of the Dirac equations.
Definition: nucleus_rmf.h:240
double energy
Energy eigenvalue.
Definition: nucleus_rmf.h:272
double omega0
Vector field at r=0.
Definition: nucleus_rmf.h:490
double rho_rhs(double sig, double ome, double rho)
Isubvector density RHS.
void meson_init()
Calculate meson and photon Green's functions gin and gout.
double rprms
Proton RMS radius (in fm)
Definition: nucleus_rmf.h:374
shell proton_shells[n_internal_levels]
The starting proton levels.
Definition: nucleus_rmf.h:548
ubvector xrhor
The isubvector field RHS.
Definition: nucleus_rmf.h:678
ubmatrix field0
Values of the fields from the last iteration.
Definition: nucleus_rmf.h:563
Spherical closed-shell nuclei with a relativistic mean-field model in the Hartree approximation...
Definition: nucleus_rmf.h:230
double etot
Total energy (in MeV)
Definition: nucleus_rmf.h:381
bool err_nonconv
If true, call the error handler if the routine does not converge or reach the desired tolerance (defa...
Definition: nucleus_rmf.h:443
static const int grid_size
The grid size.
Definition: nucleus_rmf.h:551
ubvector chdenc
Charge density.
Definition: nucleus_rmf.h:684
double A0
Coulomb field at r=0.
Definition: nucleus_rmf.h:494
double step_size
The grid step size (default 0.04)
Definition: nucleus_rmf.h:554
double r_charge_cm
Charge radius corrected by the center of mass (in fm)
Definition: nucleus_rmf.h:393
ode_rkck_gsl< ubvector, ubvector, ubvector, ode_funct11 > def_step
The default stepper.
Definition: nucleus_rmf.h:748
int nlevels
The number of levels.
Definition: nucleus_rmf.h:329
double sigma0
Scalar field at r=0.
Definition: nucleus_rmf.h:488
void meson_solve()
Solve for the meson profiles.
A shell of nucleons for nucleus_rmf.
Definition: nucleus_rmf.h:266
o2_shared_ptr< table_units<> >::type profiles
The radial profiles.
Definition: nucleus_rmf.h:529
void init_meson_density()
Initialize the meson and photon fields, the densities, etc.
ubmatrix gin
The Green's functions inside.
Definition: nucleus_rmf.h:569
ode_step< ubvector, ubvector, ubvector, ode_funct11 > * ostep
The ODE stepper.
Definition: nucleus_rmf.h:752
double dirac_tol2
Second tolerance for Dirac equations (default ).
Definition: nucleus_rmf.h:471
ubvector xrhos
The scalar field RHS.
Definition: nucleus_rmf.h:672
double match_point
Matching radius (in fm)
Definition: nucleus_rmf.h:278
double fermi_radius
The radius for which the fields are half their central value.
Definition: nucleus_rmf.h:496
std::vector< shell > levels
The levels (protons first, then neutrons)
Definition: nucleus_rmf.h:338
void meson_iter(double ic)
Calculate meson and photon fields.
Relativistic mean field theory EOS.
Definition: eos_had_rmf.h:294
int set_eos(eos_had_rmf &r)
Set the base EOS to be used.
Definition: nucleus_rmf.h:411
interp_vec< ubvector > * gi
Interpolation object.
Definition: nucleus_rmf.h:718
double rnrp
Skin thickness (in fm)
Definition: nucleus_rmf.h:360
shell neutron_shells[n_internal_levels]
The starting neutron levels.
Definition: nucleus_rmf.h:545
double meson_tol
Tolerance for meson field equations (default ).
Definition: nucleus_rmf.h:474
ubmatrix fields
The values of the fields.
Definition: nucleus_rmf.h:566
ubvector energy
Energy integrand.
Definition: nucleus_rmf.h:691
bool init_called
True if init() has been called.
Definition: nucleus_rmf.h:762
o2_shared_ptr< table_units<> >::type chden_table
The final charge densities.
Definition: nucleus_rmf.h:533
double xpform(double x, double xp, double a)
Function representing proton form factor.
int itmax
Maximum number of total iterations (default 70)
Definition: nucleus_rmf.h:453
eos_had_rmf def_rmf
The default equation of state (default NL3)
Definition: nucleus_rmf.h:404
int load_nl3(eos_had_rmf &r)
Load the default model NL3 into the given eos_had_rmf object.
void init_run(int nucleus_Z, int nucleus_N, int unocc_Z, int unocc_N)
Initialize a run.
ubvector xrhosp
The proton scalar density times radius squared.
Definition: nucleus_rmf.h:669
void set_step(ode_step< ubvector, ubvector, ubvector, ode_funct11 > &step)
Set the stepper for the Dirac differential equation.
Definition: nucleus_rmf.h:446
double rho0
Isubvector field at r=0.
Definition: nucleus_rmf.h:492
void center_mass_corr(double atot)
Compute the center of mass correction.
int odefun(double x, size_t nv, const ubvector &y, ubvector &dydx, odparms &op)
The form of the Dirac equations for the ODE stepper.
fermion n
The neutron.
Definition: nucleus_rmf.h:427
initial_guess ig
Parameters for initial guess.
Definition: nucleus_rmf.h:505
double fermi_width
The "width" of the Fermi-Dirac function.
Definition: nucleus_rmf.h:498
double eigen
Eigenvalue.
Definition: nucleus_rmf.h:242
std::vector< shell > * levp
A pointer to the current vector of levels (either levels or unocc_levels)
Definition: nucleus_rmf.h:538
double a_proton
The parameter for the charge density of the proton (default is about 4.27073)
Definition: nucleus_rmf.h:519
void gauss(double xmin, double xmax, double x, double &xi)
Perform integrations for form factor.
double sigma_rhs(double sig, double ome, double rho)
Scalar density RHS.
double omega_rhs(double sig, double ome, double rho)
Vector density RHS.
ubvector arho
Baryon density.
Definition: nucleus_rmf.h:687
int surf_index
The grid index corresponding to the nuclear surface (computed by init_run())
Definition: nucleus_rmf.h:659
int post_converge(int nucleus_Z, int nucleus_N, int unocc_Z, int unocc_N)
After convergence, make CM corrections, etc.
int meson_itmax
Maximum number of iterations for solving the meson field equations (default 10000) ...
Definition: nucleus_rmf.h:458
ubvector xrhov
The vector field RHS.
Definition: nucleus_rmf.h:675
ubvector ode_dydx
ODE derivatives.
Definition: nucleus_rmf.h:768
void field(double x, double &s, double &v, ubvector &varr)
Compute the fields for the Dirac equations.
thermo hb
thermo object for the EOS
Definition: nucleus_rmf.h:420
double isospin
Isospin ( or .
Definition: nucleus_rmf.h:274
Initial guess structure.
Definition: nucleus_rmf.h:486
ubvector chden1
Charge density.
Definition: nucleus_rmf.h:681
ubvector ode_yerr
ODE errors.
Definition: nucleus_rmf.h:771
int dirac_itmax
Maximum number of iterations for solving the Dirac equations (default 100)
Definition: nucleus_rmf.h:463
void set_verbose(int v)
Set output level.
Definition: nucleus_rmf.h:298
o2_shared_ptr< table_units<> >::type get_profiles()
Get the radial profiles.
Definition: nucleus_rmf.h:325
double rnrms
Neutron RMS radius (in fm)
Definition: nucleus_rmf.h:367
int run_nucleus(int nucleus_Z, int nucleus_N, int unocc_Z, int unocc_N)
Computes the structure of a nucleus with the specified number of levels.
ubvector ode_y
ODE functions.
Definition: nucleus_rmf.h:765
o2_shared_ptr< table_units<> >::type get_chden()
The final charge densities.
Definition: nucleus_rmf.h:329
double dirac_rk4(double x, double g1, double f1, double &funt, double eigen, double kappa, ubvector &varr)
Integrate the Dirac equations using a simple inline 4th order Runge-Kutta.
std::string state
Angular momentum-spin state .
Definition: nucleus_rmf.h:276
void energies(double xpro, double xnu, double e)
Calculate the energy profile.
double r_charge
Charge radius (in fm)
Definition: nucleus_rmf.h:387
double xrhop(double x1)
Desc.
int nodes
Number of nodes in the wave function.
Definition: nucleus_rmf.h:284
double mnuc
The nucleon mass (automatically set in init_fun())
Definition: nucleus_rmf.h:557
fermion p
The proton.
Definition: nucleus_rmf.h:434
void dirac_step(double &x, double h, double eigen, double kappa, ubvector &varr)
Take a step in the Dirac equations.
ubmatrix xrho
The densities times radius squared.
Definition: nucleus_rmf.h:666
int nuolevels
The number of unoccupied levels (equal to unocc_Z + unocc_N)
Definition: nucleus_rmf.h:341
int verbose
Control output (default 1)
Definition: nucleus_rmf.h:542
std::vector< shell > unocc_levels
The unoccupied levels (protons first, then neutrons)
Definition: nucleus_rmf.h:347
double kappa
Quantum number .
Definition: nucleus_rmf.h:244
static const int n_internal_levels
The total number of shells stored internally.
Definition: nucleus_rmf.h:252
std::function< int(double, size_t, const boost::numeric::ublas::vector< double > &, boost::numeric::ublas::vector< double > &)> ode_funct11
void iterate(int nucleus_Z, int nucleus_N, int unocc_Z, int unocc_N, int &iconverged)
Perform an iteration.
ubmatrix * fields
The meson fields.
Definition: nucleus_rmf.h:246
ubmatrix gout
The Green's functions outside.
Definition: nucleus_rmf.h:572

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).
Hosted at Get Object-oriented Scientific Computing
Lib at SourceForge.net. Fast, secure and Free Open Source software
downloads..