Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
nf_angularMomentumCoupling.cc File Reference
#include <stdlib.h>
#include <cmath>
#include "nf_specialFunctions.h"

Go to the source code of this file.

Functions

double nf_amc_log_factorial (int n)
 
double nf_amc_factorial (int n)
 
double nf_amc_wigner_3j (int j1, int j2, int j3, int j4, int j5, int j6)
 
double nf_amc_wigner_6j (int j1, int j2, int j3, int j4, int j5, int j6)
 
double nf_amc_wigner_9j (int j1, int j2, int j3, int j4, int j5, int j6, int j7, int j8, int j9)
 
double nf_amc_racah (int j1, int j2, int l2, int l1, int j3, int l3)
 
double nf_amc_clebsh_gordan (int j1, int j2, int m1, int m2, int j3)
 
double nf_amc_z_coefficient (int l1, int j1, int l2, int j2, int s, int ll)
 
double nf_amc_zbar_coefficient (int l1, int j1, int l2, int j2, int s, int ll)
 
double nf_amc_reduced_matrix_element (int lt, int st, int jt, int l0, int j0, int l1, int j1)
 

Function Documentation

◆ nf_amc_clebsh_gordan()

double nf_amc_clebsh_gordan ( int  j1,
int  j2,
int  m1,
int  m2,
int  j3 
)

Definition at line 288 of file nf_angularMomentumCoupling.cc.

288 {
289/*
290* Clebsh-Gordan coefficient
291* = <j1,j2,m1,m2|j3,m1+m2>
292* = (-)^(j1-j2+m1+m2) * std::sqrt(2*j3+1) * / j1 j2 j3 \
293* \ m1 m2 -m1-m2 /
294*
295* Note: Last value m3 is preset to m1+m2. Any other value will evaluate to 0.0.
296*/
297
298 int m3, x1, x2, x3, y1, y2, y3;
299 double cg = 0.0;
300
301 if ( j1 < 0 || j2 < 0 || j3 < 0) return( 0.0 );
302 if ( j1 + j2 + j3 > 2 * MAX_FACTORIAL ) return( INFINITY );
303
304 m3 = m1 + m2;
305
306 if ( ( x1 = ( j1 + m1 ) / 2 + 1 ) <= 0 ) return( 0.0 );
307 if ( ( x2 = ( j2 + m2 ) / 2 + 1 ) <= 0 ) return( 0.0 );
308 if ( ( x3 = ( j3 - m3 ) / 2 + 1 ) <= 0 ) return( 0.0 );
309
310 if ( ( y1 = x1 - m1 ) <= 0 ) return( 0.0 );
311 if ( ( y2 = x2 - m2 ) <= 0 ) return( 0.0 );
312 if ( ( y3 = x3 + m3 ) <= 0 ) return( 0.0 );
313
314 if ( j3 == 0 ){
315 if ( j1 == j2 ) cg = ( 1.0 / std::sqrt( (double)j1 + 1.0 ) * ( ( y1 % 2 == 0 ) ? -1:1 ) );
316 }
317 else if ( (j1 == 0 || j2 == 0 ) ){
318 if ( ( j1 + j2 ) == j3 ) cg = 1.0;
319 }
320 else {
321 if( m3 == 0 && std::abs( m1 ) <= 1 ){
322 if( m1 == 0 ) cg = cg1( x1, x2, x3 );
323 else cg = cg2( x1 + y1 - y2, x3 - 1, x1 + x2 - 2, x1 - y2, j1, j2, j3, m2 );
324 }
325 else if ( m2 == 0 && std::abs( m1 ) <=1 ){
326 cg = cg2( x1 - y2 + y3, x2 - 1, x1 + x3 - 2, x3 - y1, j1, j3, j3, m1 );
327 }
328 else if ( m1 == 0 && std::abs( m3 ) <= 1 ){
329 cg = cg2( x1, x1 - 1, x2 + x3 - 2, x2 - y3, j2, j3, j3, -m3 );
330 }
331 else cg = cg3( x1, x2, x3, y1, y2, y3 );
332 }
333
334 return( cg );
335}

Referenced by nf_amc_reduced_matrix_element(), nf_amc_wigner_3j(), nf_amc_z_coefficient(), and nf_amc_zbar_coefficient().

◆ nf_amc_factorial()

double nf_amc_factorial ( int  n)

Definition at line 96 of file nf_angularMomentumCoupling.cc.

96 {
97/*
98* returns n! for pre-computed table. INFINITY is return if n is negative or too large.
99*/
100 return G4Exp( nf_amc_log_factorial( n ) );
101}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
double nf_amc_log_factorial(int n)

◆ nf_amc_log_factorial()

double nf_amc_log_factorial ( int  n)

Definition at line 85 of file nf_angularMomentumCoupling.cc.

85 {
86/*
87* returns ln( n! ).
88*/
89 if( n > MAX_FACTORIAL ) return( INFINITY );
90 if( n < 0 ) return( INFINITY );
91 return nf_amc_log_fact[n];
92}

Referenced by nf_amc_factorial().

◆ nf_amc_racah()

double nf_amc_racah ( int  j1,
int  j2,
int  l2,
int  l1,
int  j3,
int  l3 
)

Definition at line 253 of file nf_angularMomentumCoupling.cc.

253 {
254/*
255* Racah coefficient definition in Edmonds (AR Edmonds, "Angular Momentum in Quantum Mechanics", Princeton (1980) is
256* W(j1, j2, l2, l1 ; j3, l3) = (-1)^(j1+j2+l1+l2) * { j1 j2 j3 }
257* { l1 l2 l3 }
258* The call signature of W(...) appears jumbled, but hey, that's the convention.
259*
260* This convention is exactly that used by Blatt-Biedenharn (Rev. Mod. Phys. 24, 258 (1952)) too
261*/
262
263 double sig;
264
265 sig = ( ( ( j1 + j2 + l1 + l2 ) % 4 == 0 ) ? 1.0 : -1.0 );
266 return sig * nf_amc_wigner_6j( j1, j2, j3, l1, l2, l3 );
267}
double nf_amc_wigner_6j(int j1, int j2, int j3, int j4, int j5, int j6)

Referenced by nf_amc_wigner_9j(), nf_amc_z_coefficient(), and nf_amc_zbar_coefficient().

◆ nf_amc_reduced_matrix_element()

double nf_amc_reduced_matrix_element ( int  lt,
int  st,
int  jt,
int  l0,
int  j0,
int  l1,
int  j1 
)

Definition at line 473 of file nf_angularMomentumCoupling.cc.

473 {
474/*
475* Reduced Matrix Element for Tensor Operator
476* = < l1j1 || T(YL,sigma_S)J || l0j0 >
477*
478* M.B.Johnson, L.W.Owen, G.R.Satchler
479* Phys. Rev. 142, 748 (1966)
480* Note: definition differs from JOS by the factor sqrt(2j1+1)
481*/
482 int llt;
483 double x1, x2, x3, reduced_mat, clebsh_gordan;
484
485 if ( parity( lt ) != parity( l0 ) * parity( l1 ) ) return( 0.0 );
486 if ( std::abs( l0 - l1 ) > lt || ( l0 + l1 ) < lt ) return( 0.0 );
487 if ( std::abs( ( j0 - j1 ) / 2 ) > jt || ( ( j0 + j1 ) / 2 ) < jt ) return( 0.0 );
488
489 llt = 2 * lt;
490 jt *= 2;
491 st *= 2;
492
493 if( ( clebsh_gordan = nf_amc_clebsh_gordan( j1, j0, 1, -1, jt ) ) == INFINITY ) return( INFINITY );
494
495 reduced_mat = 1.0 / std::sqrt( 4 * M_PI ) * clebsh_gordan / std::sqrt( jt + 1.0 ) /* BRB jt + 1 <= 0? */
496 * std::sqrt( ( j0 + 1.0 ) * ( j1 + 1.0 ) * ( llt + 1.0 ) )
497 * parity( ( j1 - j0 ) / 2 ) * parity( ( -l0 + l1 + lt ) / 2 ) * parity( ( j0 - 1 ) / 2 );
498
499 if( st == 2 ){
500 x1 = ( l0 - j0 / 2.0 ) * ( j0 + 1.0 );
501 x2 = ( l1 - j1 / 2.0 ) * ( j1 + 1.0 );
502 if ( jt == llt ){
503 x3 = ( lt == 0 ) ? 0 : ( x1 - x2 ) / std::sqrt( lt * ( lt + 1.0 ) );
504 }
505 else if ( jt == ( llt - st ) ){
506 x3 = ( lt == 0 ) ? 0 : -( lt + x1 + x2 ) / std::sqrt( lt * ( 2.0 * lt + 1.0 ) );
507 }
508 else if ( jt == ( llt + st ) ){
509 x3 = ( lt + 1 - x1 - x2 ) / std::sqrt( ( 2.0 * lt + 1.0 ) * ( lt + 1.0 ) );
510 }
511 else{
512 x3 = 1.0;
513 }
514 }
515 else x3 = 1.0;
516 reduced_mat *= x3;
517
518 return( reduced_mat );
519}
#define M_PI
Definition: SbMath.h:33
double nf_amc_clebsh_gordan(int j1, int j2, int m1, int m2, int j3)

◆ nf_amc_wigner_3j()

double nf_amc_wigner_3j ( int  j1,
int  j2,
int  j3,
int  j4,
int  j5,
int  j6 
)

Definition at line 105 of file nf_angularMomentumCoupling.cc.

105 {
106/*
107* Wigner's 3J symbol (similar to Clebsh-Gordan)
108* = / j1 j2 j3 \
109* \ j4 j5 j6 /
110*/
111 double cg;
112
113 if( ( j4 + j5 + j6 ) != 0 ) return( 0.0 );
114 if( ( cg = nf_amc_clebsh_gordan( j1, j2, j4, j5, j3 ) ) == 0.0 ) return ( 0.0 );
115 if( cg == INFINITY ) return( cg );
116 return( ( ( ( j1 - j2 - j6 ) % 4 == 0 ) ? 1.0 : -1.0 ) * cg / std::sqrt( j3 + 1.0 ) ); /* BRB j3 + 1 <= 0? */
117}

◆ nf_amc_wigner_6j()

double nf_amc_wigner_6j ( int  j1,
int  j2,
int  j3,
int  j4,
int  j5,
int  j6 
)

Definition at line 121 of file nf_angularMomentumCoupling.cc.

121 {
122/*
123* Wigner's 6J symbol (similar to Racah)
124* = { j1 j2 j3 }
125* { j4 j5 j6 }
126*/
127 int i, x[6];
128
129 x[0] = j1; x[1] = j2; x[2] = j3; x[3] = j4; x[4] = j5; x[5] = j6;
130 for( i = 0; i < 6; i++ ) if ( x[i] == 0 ) return( w6j0( i, x ) );
131
132 return( w6j1( x ) );
133}

Referenced by nf_amc_racah().

◆ nf_amc_wigner_9j()

double nf_amc_wigner_9j ( int  j1,
int  j2,
int  j3,
int  j4,
int  j5,
int  j6,
int  j7,
int  j8,
int  j9 
)

Definition at line 226 of file nf_angularMomentumCoupling.cc.

226 {
227/*
228* Wigner's 9J symbol
229* / j1 j2 j3 \
230* = | j4 j5 j6 |
231* \ j7 j8 j9 /
232*
233*/
234 int i, i0, i1;
235 double rac;
236
237 i0 = max3( std::abs( j1 - j9 ), std::abs( j2 - j6 ), std::abs( j4 - j8 ) );
238 i1 = min3( ( j1 + j9 ), ( j2 + j6 ), ( j4 + j8 ) );
239
240 rac = 0.0;
241 for ( i = i0; i <= i1; i += 2 ){
242 rac += nf_amc_racah( j1, j4, j9, j8, j7, i )
243 * nf_amc_racah( j2, j5, i, j4, j8, j6 )
244 * nf_amc_racah( j9, i, j3, j2, j1, j6 ) * ( i + 1 );
245 if( rac == INFINITY ) return( INFINITY );
246 }
247
248 return( ( ( (int)( ( j1 + j3 + j5 + j8 ) / 2 + j2 + j4 + j9 ) % 4 == 0 ) ? 1.0 : -1.0 ) * rac );
249}
double nf_amc_racah(int j1, int j2, int l2, int l1, int j3, int l3)

◆ nf_amc_z_coefficient()

double nf_amc_z_coefficient ( int  l1,
int  j1,
int  l2,
int  j2,
int  s,
int  ll 
)

Definition at line 437 of file nf_angularMomentumCoupling.cc.

437 {
438/*
439* Biedenharn's Z-coefficient coefficient
440* = Z(l1 j1 l2 j2 | S L )
441*/
442 double z, clebsh_gordan = nf_amc_clebsh_gordan( l1, l2, 0, 0, ll ), racah = nf_amc_racah( l1, j1, l2, j2, s, ll );
443
444 if( ( clebsh_gordan == INFINITY ) || ( racah == INFINITY ) ) return( INFINITY );
445 z = ( ( ( -l1 + l2 + ll ) % 8 == 0 ) ? 1.0 : -1.0 )
446 * std::sqrt( l1 + 1.0 ) * std::sqrt( l2 + 1.0 ) * std::sqrt( j1 + 1.0 ) * std::sqrt( j2 + 1.0 ) * clebsh_gordan * racah;
447
448 return( z );
449}

◆ nf_amc_zbar_coefficient()

double nf_amc_zbar_coefficient ( int  l1,
int  j1,
int  l2,
int  j2,
int  s,
int  ll 
)

Definition at line 453 of file nf_angularMomentumCoupling.cc.

453 {
454/*
455* Lane & Thomas's Zbar-coefficient coefficient
456* = Zbar(l1 j1 l2 j2 | S L )
457* = (-i)^( -l1 + l2 + ll ) * Z(l1 j1 l2 j2 | S L )
458*
459* Lane & Thomas Rev. Mod. Phys. 30, 257-353 (1958).
460* Note, Lane & Thomas define this because they did not like the different phase convention in Blatt & Biedenharn's Z coefficient. They changed it to get better time-reversal behavior.
461* Froehner uses Lane & Thomas convention as does T. Kawano.
462*/
463 double zbar, clebsh_gordan = nf_amc_clebsh_gordan( l1, l2, 0, 0, ll ), racah = nf_amc_racah( l1, j1, l2, j2, s, ll );
464
465 if( ( clebsh_gordan == INFINITY ) || ( racah == INFINITY ) ) return( INFINITY );
466 zbar = std::sqrt( l1 + 1.0 ) * std::sqrt( l2 + 1.0 ) * std::sqrt( j1 + 1.0 ) * std::sqrt( j2 + 1.0 ) * clebsh_gordan * racah;
467
468 return( zbar );
469}