DynaRho

deMon-DynaRho v. 3.3 Manual Last update: 11 April 2002

If you experience difficulties with this page, please contact me at Mark.Casida@ujf-grenoble.fr.

Manual for deMon-DynaRho version 3.3 (xxx release)

Neptune and dolphin

(Gian Lorenzo Bernini, Italian, Rome, about 1620-1680)

Cite this program as: deMon-DynaRho version 3.3, M.E. Casida, C. Jamorski, E. Fadda, J. Guan, S. Hamel, and D.R. Salahub, University of Montreal and Université Joseph-Fourier, copyright 2003.

Table of Contents

[1.] OVERVIEW [2.] USERS' GUIDE [2.1] INPUT [2.2] OUTPUT [3.] PROGRAMMER'S GUIDE [3.1] LANGUAGE AND SUPPORTED PLATFORMS [3.2] PROGRAM STRUCTURE [3.3] UTILITY ROUTINES [3.4] DYNAMICAL MEMORY ALLOCATION [3.5] FILES [4.] BIBLIOGRAPHY

[1.] OVERVIEW

This part of this manual is now in

in order to facilitate typing equations and reference handling.

[2.] USER'S GUIDE

DynaRho version 3 is under rapid development. For this reason, the input and output is changing rapidly. Thus far, the program is restricted to closed-shell systems only, where the spin alpha and spin beta orbitals are identical. Version 3.2 now calculates its own grid and integrals, making it relatively stand alone. As such it can work with many different versions of deMon-KS. Schematically, the file transfer between deMon-KS and deMon-DynaRho is now,

+---------------+ BASIS ---------------------------------> | | | +----------+ | | +------> | | -> ALCHEM.inp -> | deMon-DynaRho | -> Dyna.out | deMon-KS | | | SCF.inp -> | | Dyna.inp ---> | | +----------+ +---------------+

Compatible versions of deMon-KS include:

  • deMon-KS3p4 --- Note that the "ALCHEM" keyword must be used. "NOFSYM" must be used, as well as all six cartessian d-orbitals (the default).
  • deMon-KS4p0 --- Note that the "ALCHEM" and "ORBI 6D" keywords must be used. The use of "FSYM" is OK.

[2.1] INPUT

The input file is divided into two parts, a preliminary and a main part.

Preliminary input

The preliminary input consists of the

  1. atomic number
  2. auxiliary basis name
  3. orbital basis name

for each atom, in the order that the atoms appear in the SCF output, followed by the name of the grid to be used in the response calculation. These parameters are normally the same as in the SCF step.

Example:

9 'A-FLUORINE (4,4;4,4)' 'O-FLUORINE (621/41/1*)' 1 'A-HYDROGEN (3,1;3,1)' 'O-HYDROGEN (41/1*)' 'GRID FINE RANDOM'

Main input

The preliminary input is followed by the main input, described here one by one. The minimum amount of input is indicated using capital letters, while small letters indicate the rest of the keyword. When one or more options are possible for a given keyword, they are indicated as choices within parentheses.

  1. MODULE RESPONSE
    1. This is always the first line of the input file.
  2. EXCItation (FULL, DAVIdson nwant rcutoff limit maxiter)
    1. The EXCItation keyword requests the calculation of the excitation spectrum. The default is not to calculate the spectrum.
    2. The option FULL means to set up and diagonalize the full response theory problem. This is possible for only the smallest molecules.
    3. The option DAVIdson allows larger molecules to be treated by using the block Davidson algorithm for finding the lowest eigenvalues and eigenvectors of a large matrix (without actually creating and storing the large matrix). The DAVIdson keyword must be followed by the four numbers "nwant rcutoff limit maxiter," where
      • nwant is (half) the number of guesses which the block Davidson method will try to converge.
      • rcutoff is the convergence criterion (in eV) for the excitation energies.
      • limit is the dimension of the largest matrices which should be diagonalized.
      • maxiter is the maximum number of Davidson iterations to be allowed.
    4. BE CAREFUL: There are no defaults for these four numbers!
  3. INTEgrals (MO, AO) (ACTIve minorb maxorb)
    1. Form the coupling matrix or its action on a vector by either an AO-based or MO-based algorithm. The default is MO.
    2. The ACTIve keyword is optional and restricts calculations to a subspace of active molecular orbitals. Here
      • minorb is the number of the smallest active orbital.
      • maxorb is the number of the largest active orbital.
  4. POLArizability (STATic, DYNAmic)
    1. The POLArizability key word requests the calculation of the dipole polarizability by matrix inversion. The default is not to do the calculation. Either STATic or DYNAmic polarizabilities may be calculated.
  5. PREScreening THR1 THR2
    1. This sets a threshold of 10**(-THR1) for local prescreening in the calculation of the overlap and exchange-correlation second-derivative matrices in the auxiliary basis representation, and a threshold of 10**(-THR2) for global prescreening used in the construction of the exchange-correlation B matrices. The default is THR1=8 and THR2=8. For more information and test results click here.
  6. SCF (LDA, HF, OEP, KLI) (FIRST, (MAXSCF,ENCONV,DMIX)) The SCF keyword specifies deMon-DynaRho should carry out the self-consistent field calculation which preceeds the response theory calculation, instead of just taking the result directly from deMon-KS. The choices of SCF calculation in DynaRho are:
    • LDA --- local density approximation
    • HF --- Hartree-Fock
    • OEP --- exchange-only optimized effective potential
    • KLI --- pseudo Krieger-Li-Iafrate (i.e. the average denominator approximation is made in the OEP calculation)
    1. All of these are done in the resolution-of-the-identity (RI) approximation. These keywords need to be followed by either:
      • FIRST --- first-order calculation only (i.e. single iteration)
      • (MAXSCF,ENCONV,DMIX) --- iterate with a mixing of DMIX until an total energy change of ENCONV for eight consecutive iterations for a maximum of MAXSCF iterations
  7. LEVEl (IPA, RPA, XA, LDA)
    1. The LEVEl keyword specifies the level of calculation of the TD-DFRT coupling matrix. Here, IPA denotes the independent particle approximation, RPA denotes the density-functional theory random phase approximation, XA denotes the time-dependent X alpha method with alpha=2/3, while LDA denotes the full time-dependent local density approximation (TDLDA).
  8. TWOLevel (NMR)
    1. For closed-shell TDLDAxc calculations, this calculates the singlet and triplet TDLDA and DeltaSCF energies in the twolevel model described in Ref. [C96]. If the NMR option is used, then an NMR input file will be created for a modified version of the chemical shift part of deMon-NMR.
  9. SPINblocking (SINGlet,TRIPlet,BOTH,NONE)
    1. Explicitly block the TDDFT Omega matrix according to spin-coupling.
      • SINGlet Calculate singlet-coupled results only.
      • TRIPlet Calculate triplet-coupled results only.
      • BOTH The default. Calculate both singlet- and triplet-coupled results.
      • NONE Do not block according to spin-coupling.
  10. PROPERTY DIPOLE
    1. This line is mandatory.
  11. END
    1. This is always the last line.

PRE> Omega F_I = omega_I^2 F_I , is solved. Here omega_I is the excitation energy for the Ith transition and the F_I is the coefficient vector for the transition density. This vector is normalized to one,

|F_I| = 1 ,

and is related to the transition densities via the equations,

f_I^x = |x^+ S^(-1/2) F_I |^2 , f_I^y = |y^+ S^(-1/2) F_I |^2 , f_I^z = |z^+ S^(-1/2) F_I |^2 , f_I = (2/3)(f_I^x + f_I^y + f_I^z ) ,

where

S^(-1/2)_{ia,jb} = delta_{i,j} delta_{a,b} * (epsilon_a - epsilon_i )^{1/2}

is the diagonal matrix corresponding to the squareroot of the IPA excitation energies.

Explanation of the table :

+-------------------------------+ | A negative excitation energy | | denotes an imaginary energy. | +-------------------------------+ --------------------------------------------------------------------- ! TRANSITION ! TRANS ENERGY ! OSCIL STRENGTH ! fi/wi2 ! --------------------------------------------------------------------- ! 1 ! 0.429423 ! 0.000000 ! 0.000000 ! ! 2 ! 0.755784 ! 0.000000 ! 0.000000 ! ! 3 ! 0.850347 ! 0.000000 ! 0.000000 ! ! 4 ! 1.252958 ! 0.000000 ! 0.000000 ! ! 5 ! 1.408248 ! 0.000000 ! 0.000000 ! ! 6 ! 1.590388 ! 0.000000 ! 0.000000 ! ! 7 ! 1.824482 ! 0.000000 ! 0.000000 ! ! 8 ! 1.933882 ! 0.000101 ! 0.020003 ! ! 9 ! 1.962512 ! 1.425585 ! 274.061273 ! ! 10 ! 2.083548 ! 0.000000 ! 0.000000 !

Column 1 gives the order of the transitions. Column 2 gives omega_I in eV. Since omega_I^2 may sometimes be negative (when this occurs it may be related to instability of the spin restricted solution), omega_I may turn out to be imaginary. For convenience of printing,

sign(omega_I^2) sqrt(abs(omega_I^2))

is actually printed instead of omega_I . Column 3 is the oscillator strength (unitless). Column 4 gives f_I/omega_I^2 in atomic units and is useful because the mean polarizability is given by the sum over states expression,

alpha = Sum_I f_I/omega_I^2 .

The Cauchy coefficients are defined by

S_J = Sum_I (omega_I)^J f_I .

The Cauchy coefficients are defined by

S_J = Sum_I (omega_I)^J f_I .

Some Cauchy coefficients have a definite interpretation,

S_0 in the limit of a complete basis set = the number of electrons (TRK sum rule) , S_{-2} = mean static polarizability .

They are printed in a.u. in a table :

*-------------------------------------------* | Oscillator Strength Sums | *-------------------------------------------* (Frozen-orbital-approximation sums have been refined where possible by Davidson iterations.) S_2 0.742141043E-02 S_1 0.102885272 S_0 (TRK sum) 1.42636826 S_-2 (mean polarizability) 274.166970 S_-4 52701.5760 S_-6 10130909.2 S_-8 0.194752899E+10 S_-10 0.374391895E+12 S_-12 0.719736504E+14 S_-14 0.138364171E+17

Note that these are not accurate for frozen orbital or Davidson calculations, and are only presented as (usually bad) first estimate of the true Cauchy coefficients.

The excited state symmetry may be determined from the transition density coefficients which are printed in a table :

!------------------------------------------------------------! Transition 1 Energy 3.043270 !------------------------------------------------------------! ! 8 --> 9 ! -0.247666 ! 2 B1 --> 2 B2 ! !------------------------------------------------------------! ! 8 --> 9 ! 0.247666 ! 2 B1 --> 2 B2 ! !------------------------------------------------------------! Transition 2 Energy 3.647364 !------------------------------------------------------------! ! 8 --> 9 ! -0.247651 ! 2 B1 --> 2 B2 ! !------------------------------------------------------------! ! 8 --> 9 ! -0.247651 ! 2 B1 --> 2 B2 ! !------------------------------------------------------------! Transition 3 Energy 5.757730 !------------------------------------------------------------! ! 8 --> 10 ! -0.326817 ! 2 B1 --> 6 A1 ! ! 8 --> 12 ! 0.013372 ! 2 B1 --> 7 A1 ! ! 8 --> 16 ! -0.013664 ! 2 B1 --> 9 A1 ! ! 8 --> 25 ! 0.010584 ! 2 B1 --> 13 A1 ! !------------------------------------------------------------! ! 8 --> 10 ! 0.326817 ! 2 B1 --> 6 A1 ! ! 8 --> 12 ! -0.013372 ! 2 B1 --> 7 A1 ! ! 8 --> 16 ! 0.013664 ! 2 B1 --> 9 A1 ! ! 8 --> 25 ! -0.010584 ! 2 B1 --> 13 A1 ! !------------------------------------------------------------! Transition 4 Energy 5.832088 !------------------------------------------------------------! ! 8 --> 10 ! 0.327586 ! 2 B1 --> 6 A1 ! !------------------------------------------------------------! ! 8 --> 10 ! 0.327586 ! 2 B1 --> 6 A1 ! !------------------------------------------------------------!

The automatic orbital assignments are only available when DynaRho is run with KS4p0, otherwise the orbitals have to be assigned by hand. The transition energy is in eV. The coefficients printed are actually those components of

S^{-1/2} F_I

which exceed 0.01 a.u. Direct comparison with CI coefficients requires division of these coefficients by omega_I (in a.u.) (see Ref. [C95] pg. 402). The coefficients have also been divided into spin alpha components (above the dashed line) and spin beta components (below). This makes it easy to assign singlet and triplet states in a closed-shell molecule since the spin alpha and beta coefficients have the opposite sign for triplet states (as in Transition 1) and the same sign for singlet states (as in Transition 2).

Tamm-Dancoff Approximation

As we have seen the TDA equation is

A XI = omegaI XI .

The eigenvalue, omega_I, is still the excitation energy for the Ith transition, but oscillator strengths are now calculated differently. The eigenvectors are normalized so that

|XI| = 1 ,

and are now related to the transition densities via the equations,

fIx = omegaI |x+ XI |2 , fIy = omegaI |y+ XI |2 , fIz = omegaI |z+ XI |2 , fI = (2/3)(fIx + fIy + fIz ) .

Note how the S matrix (whose elements are basically differences of orbital energies) is completely absent and has been replaced with the excitation energy, omegaI. This looks like a big change, but the effect on oscillator strengths can be surprisingly small.

[3.] PROGRAMMER'S GUIDE

Any programmer who wants to work with DynaRho is seriously encouraged to begin by mastering the equations presented in Refs. [C95,C96].

[3.1] LANGUAGE AND SUPPORTED PLATFORMS

DynaRho is written in standard FORTRAN77 with the following exceptions:

  1. some calls to system routines to get time and date information,
  2. C language syntax include statements, and
  3. two C routines to do dynamical memory allocation.

So far DynaRho has been successfully ported to the following platforms:

  • SGI R4000
  • SGI R8000: use the -Dbit64 compile option
  • Linux/pgf77: use the -Dpgf77 option. See the porting notes.

[3.2] MAIN PROGRAM STRUCTURE

MAIN ..HEADER (print header) ..INICST (initialize constants) ..INPUT (read input files) ....INIPROPS (read input file created by DynaSCF) ......INBASIS (read basis functions) ........OUTBAS (print table of basis functions) ......GRID1 (memory allocation for GRID) ........GRID (create grid) ..........DIST (function to calculate distance between atoms) ..........ARANSET (sets JRAND=IDUMMY) ..........POINTS (chooses angular points and corrects weight factor) ............ARANF (function) ............ANG50 ............ANG110 ............ANG194 ..........COMPGR (estimate importance of grid point for basis) ......SCD (driver for calculating inverse CD overlap matrix) ........TAG (initialization of arrays for incomplete gamma function) ........CDSCFSS (calculate s,s-type elements of CD overlap matrix) ..........AUGG (calculate incomplete gamma function) ........CDSCFXS (calculate spd,s-type elements of CD overlap matrix) ..........AUGG (calculate incomplete gamma function) ........CDSCFXX (calculate spd,spd-type elementsof CD overlap matrix) ..........AUGG (calculate incomplete gamma function) ........INVERS (invert the CD overlap matrix) ......FOCKSS1 (memory allocation for FOCKSS) ........FOCKSS (calculate ss-type 3-center ERIs) .........AUGG (calculate incomplete gamma function) ......FOCKPS1 (memory allocation for FOCKPS) ........FOCKPS (calculate ps-type 3-center ERIs) ..........AUGG (calculate incomplete gamma function) ......FOCKPP1 (memory allocation for FOCKPP) ........FOCKPP (calculate pp-type 3-center ERIs) ..........AUGG (calculate incomplete gamma function) ......FOCKDS1 (memory allocation for FOCKDS) ........FOCKDS (calculate ds-type 3-center ERIs) ..........AUGG (calculate incomplete gamma function) ......FOCKDP1 (memory allocation for FOCKDP) ........FOCKDP (calculate dp-type 3-center ERIs) ..........AUGG (calculate incomplete gamma function) ......FOCKDD1 (memory allocation for FOCKDD) ........FOCKDD (calculate dd-type 3-center ERIs) ..........FOCKDDA (coulomb-type dd 3-center ERIs) ...........AUGG (calculate incomplete gamma function) ......IOTHREE1 (memory allocation for IOTHREE) ........IOTHREE (read in, sort, and write out 3-center integrals) ....LECTURE (read input file created by user) ......TEXT ......DIGITS ..FORM1 (allocation of memorhy for FORM) ....FORM (transform 3-center ERIs over AO basis to MO basis) ......FILMUNU ......AOTOMO ..RIHF1 (allocation of memory for RIHF) ....RIHF (resolution-of-the identity Hartree-Fock) ......HMATLL (construct orbital overlap matrix, S, and core Hamiltonian, HMAT) ........HMATSS (s,s matrix elements) ........HMATPS (p,s matrix elements) ........HMATPP (p,p matrix elements) ........HMATDS (d,s matrix elements) ........HMATDP (d,p matrix elements) ........HMATDD (d,d matrix elements) ......GETJ1 (allocation of memory for GETJ) ........GETJ (construct coulomb potential matrix) ......OVERLL (construct orbital overlap matrix, S) ........OVERSS (s,s matrix elements) ........OVERPS (p,s matrix elements) ........OVERPP (p,p matrix elements) ........OVERDS (d,s matrix elements) ........OVERDP (d,p matrix elements) ........OVERXX (d,d matrix elements) ......BLDVXC1 (allocation of memory for BLDVXC) ........BLDVXC (constrct exchange-correlation potential) ......GETK1 (allocation of memory for GETK) ........GETK (construct HF exchange matrix) ......LOWDIN1 (allocation of memory for LOWDIN) ........LOWDIN (Solve F C = S C E ) ..........DIAG1 (allocation of memory for DIAG) ............DIAG (solve A U = U D ) ..OEP1 (allocation of memory for OEP1) ....OEP (resolution-of-the-identity optimized effective potential) ......HMATLL (construct orbital overlap matrix, S, and core Hamiltonian, HMAT) ........HMATSS (s,s matrix elements) ........HMATPS (p,s matrix elements) ........HMATPP (p,p matrix elements) ........HMATDS (d,s matrix elements) ........HMATDP (d,p matrix elements) ........HMATDD (d,d matrix elements) ......GETJ1 (allocation of memory for GETJ) ........GETJ (construct coulomb potential matrix) ......OVERLL (construct orbital overlap matrix, S) ........OVERSS (s,s matrix elements) ........OVERPS (p,s matrix elements) ........OVERPP (p,p matrix elements) ........OVERDS (d,s matrix elements) ........OVERDP (d,p matrix elements) ........OVERXX (d,d matrix elements) ......BLDVXC1 (allocation of memory for BLDVXC) ........BLDVXC (constrct exchange-correlation potential) ......GETK1 (allocation of memory for GETK) ........GETK (construct HF exchange matrix) ......GETLAMBDA1 (allocate memory for GETLAMBDA) ........GETLAMBDA (response to exchange operator) ..........FORM1 (allocation of memorhy for FORM) ............FORM (transform 3-center ERIs over AO basis to MO basis) ..............FILMUNU ..............AOTOMO ......GETX1 (allocate memory for GETX) ........GETX (get chi matrix) ..........SVDCMP (singular value decomposition from numerical recipies) ............PYTHAG ......GETB1 (allocate memory for GETB) ........GETB (solve OEP equation to obtain fitting coefficients) ......GETVXC1 (allocate memory for GETVXC) ........GETVXC (construct Vxc matrix from fitting coefficients) ......LOWDIN1 (allocation of memory for LOWDIN) ........LOWDIN (Solve F C = S C E ) ..........DIAG1 (allocation of memory for DIAG) ............DIAG (solve A U = U D ) ..KLI1 (allocation of memory for KLI) ....KLI ......HMATLL (construct orbital overlap matrix, S, and core Hamiltonian, HMAT) ........HMATSS (s,s matrix elements) ........HMATPS (p,s matrix elements) ........HMATPP (p,p matrix elements) ........HMATDS (d,s matrix elements) ........HMATDP (d,p matrix elements) ........HMATDD (d,d matrix elements) ......GETJ1 (allocation of memory for GETJ) ........GETJ (construct coulomb potential matrix) ......OVERLL (construct orbital overlap matrix, S) ........OVERSS (s,s matrix elements) ........OVERPS (p,s matrix elements) ........OVERPP (p,p matrix elements) ........OVERDS (d,s matrix elements) ........OVERDP (d,p matrix elements) ........OVERXX (d,d matrix elements) ......BLDVXC1 (allocation of memory for BLDVXC) ........BLDVXC (constrct exchange-correlation potential) ......GETK1 (allocation of memory for GETK) ........GETK (construct HF exchange matrix) ......LOWDIN1 (allocation of memory for LOWDIN) ........LOWDIN (Solve F C = S C E ) ..........DIAG1 (allocation of memory for DIAG) ............DIAG (solve A U = U D ) ..ESTIMATE (estimate memory needed) BADLY NEEDS UPDATING ..DRIVER (do calculations) ....PRINTMO (print MOs and density matrix, if desired) ....DENSTY1 (allocate memory for xc work) ......DENSTY (calculate and store the B-matrices) ........MAKEKNL (make kij and lij lookup arrays) ....GETCHI (....) ....POLAR1 (allocate memory for polar) ......POLAR (calculate polarizability) ........ACTK1 (act K on a vector) ....EXCI1 (allocate memory for exci) ......EXCI (calculate excitation spectrum) ........PRECOND1 (allocate memory for precond) ..........PRECOND (perform full or partial diagonalization) ............ACTOM1 (allocate memory for actOm) ..............ACTOM (act Omega on a vector) ................ACTK1 (act K on a vector, AO-based algorithm) ................ACTKOM1 (act K on a vector, MO-based algorithm) ............DIAG1 (Householder matrix diagonalization) ........DAVIDSON1 (allocate memory for davidson) ..........DAVIDSON (davidson diagonalization) ............ACTOM1 (allocate memory for actOm) ..............ACTOM (act Omega on a vector) ................ACTK1 (act K on a vector, AO-based algorithm) ................ACTKOM1 (act K on a vector, MO-based algorithm) ..........DIAG1 (Householder matrix diagonalization) ..........TTXYZ (Calculate matrices for X, Y, and Z) ..TWOLEVEL1 (allocate memory for TWOLEVEL) ....TWOLEVEL (TDLDAxc and DeltaSCF energies in two level model)

[3.3] UTILITY ROUTINES

ACTK1 (allocate memory for actK) ..ACTK (act K on a vector, AO-based algorithm) ....FILMUNU (fill munu lookup array) ....MOTOAO (convert vector from MO to AO rep) ....AOTOMO (convert vector from AO to MO rep) ACTKOM1 (allocate memory for actkom) ..ACTKOM (act K on a vector, MO-based algorithm) ....FILMUNU (fill munu lookup array) ....MOTOAO (convert vector from MO to AO rep) ....AOTOMO (convert vector from AO to MO rep) CLOCK (for timing) DIAG1 (allocate memory for DIAG) ..DIAG (Householder diagonalization) ....TRED2 ....TQL2 FLUSH (system call to assure printing if crash) FRAME (for printing titles) HEADER (print program header)

[3.4] DYNAMICAL MEMORY ALLOCATION

Dynamical memory allocation is accomplished in DynaRho by mixing C and FORTRAN. Such mixing is nonstandard, but apparently universally possible on UNIX systems. As is well known, dynamic memory allocation is a method to increase the amount of memory reserved for a program. While, in our experience, freeing memory does not result in a reduction of amount of memory reserved for the FORTRAN program, it does free up memory to be reused by the FORTRAN program, hence allowing adequate control of the size of the FORTRAN job. Two C functions (``fmalloc'' and ``ffree'') are used:

/* ====================== */ /* Last modified: 26/5-98 */ /* ====================== */ #include #include /* -----------------------------------------------------*/ int* fmalloc_( int*size, int*memtot) /* NOTE: The integer passed by FORTRAN is received as a */ /* pointer by C */ { /* This routine allocates memory */ int*ptr = NULL; ptr = (int *)malloc((size_t)*size); if(!ptr) { printf("Attempt to allocate %d bytes failed\n",*memtot); fflush(stdout); exit(-1); }; return ptr; } /* --------------------------------------------------------*/ void ffree_( int* *ptr ) { /* this routine de-allocates memory */ free( (void*)*ptr); return; } /* EOF */

The C functions are called from the memory allocation FORTRAN subroutines (those with names ending in the digit ``1''). For example:

c ======================================================== subroutine actk1(ispin1,Fout,ispin2,Fin) c ======================================================== c Subroutine to dynamically allocate memory for ACTK c Written: 7/12-96 c Last modified: 5/3-97 c -------------------------------------------------------- implicit real*8(a-h,o-z) c #ifdef bit64 integer*8 fmalloc integer*8 munu,faoao,faoao2,faomo,tci integer*8 Faux1,Faux2,B #else integer*4 fmalloc integer*4 munu,faoao,faoao2,faomo,tci integer*4 Faux1,Faux2,B #endif c #include "param.h" #include "common.h" c NAELEC=MAX( NLA*(NCNTRT-NLA) , NLB*(NCNTRT-NLB) ) nlmax=max(nla,nlb) c -------------------------------------------- c Calculate pointers c -------------------------------------------- mhold=memtot c nbytes=4*ncntrt*(ncntrt+1)/2 memtot=memtot+nbytes munu=fmalloc(nbytes,memtot) c nbytes=8*ncntrt*(ncntrt+1)/2 memtot=memtot+nbytes faoao=fmalloc(nbytes,memtot) c nbytes=8*ncntrt*(ncntrt+1)/2 memtot=memtot+nbytes faoao2=fmalloc(nbytes,memtot) c nbytes=8*nlmax*ncntrt memtot=memtot+nbytes faomo=fmalloc(nbytes,memtot) c npair = ncntrt*(ncntrt+1)/2 nbytes=8*nbxc*npair memtot=memtot+nbytes tci=fmalloc(nbytes,memtot) c nbytes=8*nbxc memtot=memtot+nbytes Faux1=fmalloc(nbytes,memtot) c nbytes=8*nbxc memtot=memtot+nbytes Faux2=fmalloc(nbytes,memtot) c nbytes=8*nbxc memtot=memtot+nbytes B=fmalloc(nbytes,memtot) c if(memmax.le.memtot)memmax=memtot c ------------------------------------------------------------------- call actk(ispin1,Fout,ispin2,Fin, & naelec,nlmax,npair, & %val(munu),%val(faoao),%val(faoao2), & %val(faomo),%val(tci),%val(Faux1), & %val(Faux2),%val(B)) c -------------------------------------------------------------------- c Free up the memory c -------------------------------------------------------------------- call ffree(munu) memtot=memtot-4*ncntrt*(ncntrt+1)/2 c call ffree(faoao) memtot=memtot-8*ncntrt*(ncntrt+1)/2 c call ffree(faoao2) memtot=memtot-8*ncntrt*(ncntrt+1)/2 c call ffree(faomo) memtot=memtot-8*nla*ncntrt c call ffree(tci) memtot=memtot-8*nbxc*npair c call ffree(Faux1) memtot=memtot-8*nbxc c call ffree(Faux2) memtot=memtot-8*nbxc c call ffree(B) memtot=memtot-8*nbxc c if(mhold.ne.memtot)then write(6,*)'accounting error in actk1' call flush(6) endif c --------------------------------------------------------------------- return end

In using these routines, it is important to realize that the C functions are passing and receiving pointers while the FORTRAN routines are passing and receiving integers, hence interconversion is going on. The pointers are integers whose length depends upon the word length of the compilation. This is controlled at compile time with the use of the C preprocessor (cpp) and the ``#ifdef bit64'' command. In particular, for a 64 bit compilation, the compile option "-Dbit64" must be used, otherwise pointer lengths appropriate for 32 bit compilation will be assumed.

[3.5] FILES

[4.] BIBLIOGRAPHY

This part of this manual is now in

in order to facilitate typing equations and reference handling.