CONGEN
5.0
Configuration generation for SCATCI
|
Projection on spin states. More...
Functions/Subroutines | |
subroutine, private | cntrct (nelt, no, ndo, cdo, thres) |
Throw away determinants with negligible contribution. More... | |
subroutine, private | dophz (nftw, nocsf, nelt, ndtrf, nconf, indo, ndo, lenndo, icdo, cdo, lencdo, iphz, leniphz, iphz0, leniphz0, nctarg, nctgt, notgt, mrkorb, mdegen, ntgsym, mcont, symtyp, npflg) |
Phase factor. More... | |
subroutine, private | dophz0 (nftw, nocsf, nelt, ndtrf, nconf, indo, ndo, lenndo, icdo, cdo, lencdo, iphz, npflg) |
Phase factor. More... | |
integer function, private | iphase (nconf, nelt) |
Sequence phase factor. More... | |
subroutine, private | mkorbs (nob, nsym, mn, mg, mm, ms, norb, nsrb_in, map, mpos, iposit, nobl, nob0l, symtyp) |
Compute the orbital table. More... | |
subroutine, private | pkwf (nod, ieltp, cdo, mdo, idopl, mdop, idcpl, mdcp, nftw, ndo, ndto, len_ndto, ithis_csf) |
Pack wave function. More... | |
subroutine, private | pmkorbs (nob, nobe, nsym, mn, mg, mm, ms, nsrb, norb, nsrbd, map, mpos, iposit, symtyp) |
? More... | |
subroutine, private | popnwf (nsrb, nsrbs, nelt, ndtrf, mopmx, mdop, mdcp, mop, mdc, mdo, ndta, nod, nda, idop, idcp, ieltp, flip, nalm) |
Get open-shell part of determinants. More... | |
subroutine, private | prjct (nelt, mxss, nodi, ndo, cdi, nodo, cdo, maxcdo, mgvn, iss, isd, thres, r, ndtr, mm, ms, maxndo, symtyp, nsrb) |
Apply Lowdin projection operator. More... | |
subroutine, public | projec (sname, megul, symtyp, mgvn, s, sz, r, pin, nocsf, byproj, idiag, npflg, thres, nelt, nsym, nob, ndtrf, nftw, iposit, nob0, nob1, nob01, iscat, ntgsym, notgt, nctgt, mcont, gucont, mrkorb, mdegen, mflag, nobe, nobp, nobv, maxtgsym) |
Project the wave function. More... | |
subroutine, private | ptpwf (nftw, nocsf, nelt, ndtrf, nodi, indi, icdi, ndi, cdi) |
Print the CSFs. More... | |
integer function, private | qsort (n, a) |
Sort integer array. More... | |
subroutine, private | rdwf (nft, k1, nodi, k2, cdi, k3, ndi) |
Read CSF. More... | |
subroutine, private | rdwf_getsize (iunit, num_csfs, num_dets, len_dets) |
Reads the size of the wavefunction. More... | |
subroutine, private | rfltn (nelt, nodi, ndi, cdi, r, mxnd, ndmxp, thres, nodo, cdo, ndtr, mm, bst) |
Add mirror-reflected spin-orbitals. More... | |
real(kind=wp) function, private | snrm2 (n, array, istep) |
Real strided vector norm. More... | |
subroutine, private | stmrg (nelt, maxcdo, maxndo, ndo, cdo, nodo, ndi, cdi, bst) |
Add determinant to a list of determinants. More... | |
subroutine, private | wfgntr (mgvn, iss, isd, thres, r, symtyp, nelt, nsyml, nob, nobl, nob0l, nobe, norb, nsrb, mn, mg, mm, ms, iposit, map, mpos, nocsf, ndtrf, nodi, ndi, cdi, indil, icdil, maxndi, maxcdi, nodo, ndo, cdo, indo, icdo, maxndo, maxcdo, lenndo, lencdo, npflg, byproj, nftw, nalm) |
Form spin eigenstates. More... | |
subroutine, private | wrnfto (sname, mgvn, s, sz, r, pin, norb, nsrb, nocsf, nelt, idiag, nsym, symtyp, nob, ndtrf, nodo, m, icdo, indo, ndo, lndi, cdo, lcdi, nfto, nobl, nx, npflg, thres, iposit, nob0, nob0l, nctarg, ntgsym, notgt, nctgt, mcont, gucont, iphz, nobe, nobp, nobv, maxtgsym) |
WRNFTO - WRite wavefunction data to unit NFTO. More... | |
subroutine, private | wrwf (nft, n1, nodo, n2, cdo, n3, ndo) |
Write wave functions using SPEEDY format. More... | |
Routines in this module, most prominently the central subroutine projec, read back all CSFs generated in the previous part of CONGEN execution and recombine the determinants to satisfy spin composition rules. While doing that, several new combinations of spin-orbitals (determinants) may be added if there are any open shells where the electrons spins can be freely permuted. Also, the routines apply a threshold for final selection of contributing determinants, so the output can be even smaller than the input (though this mostly signalizes some error in setup).
|
private |
Scans the store of determinants, discard such whose contribution (multiplication factor) is below given tolerance, and bubbles out the vacated intervals from the storage arrays.
Definition at line 66 of file congen_projection.f90.
Referenced by prjct(), and rfltn().
|
private |
Compute phase factor implied by placing continuum spin-orbital after all target spin-orbitals
Definition at line 107 of file congen_projection.f90.
References iphase().
Referenced by projec().
|
private |
Compute phase factor for target CSFs - given by reordering spin-orbitals in ascending order.
Definition at line 232 of file congen_projection.f90.
References iphase().
Referenced by projec().
|
private |
Compute phase factor (if any) due to out of sequence ordering of spin-orbitals in CSF stored in nconf.
Definition at line 283 of file congen_projection.f90.
Referenced by dophz(), dophz0(), and projec().
|
private |
Computes the orbital table which is then used in the projection step. This is called from projec().
Input data: ISYMTYP Switch for C-inf-v (=0 or 1) / Abelian point group (=2 NOB Number of orbitals per symmetry NSYM Number of symmetries in the orbital set NPFLG Flag controlling printing of computed orbital table Output data: MN Orbital number associated with each spin-orbital MG G/U designation for each spin-orbital (C-inf-v only) Actually this is always zero because C-inf-v does not distinguish between g/u. It exists because original version of Alchemy tried to use it for D-inf-h too; all CI evauation is doen in C-inf-v now because CONGE converts D-inf-h to C-inf-v data. MM Symmetry quantum number associated with each spin-orb MS Spin function ( alpha or beta ) associated with each spin orbital Notes: The orbital table establishes orbital and quantum number data for each spin orbital in the set. e.g. C-inf-v symmetry with NSYM=2, NOB=3,1, yields ten spin orbitals which are designated as follows by this routine: Spin orb. MN MG MM MS Comments 1 1 0 0 0 1 sigma spin up 2 1 0 0 1 1 sigma spin down 3 2 0 0 0 2 sigma spin up 4 2 0 0 1 2 sigma spin down 5 3 0 0 0 3 sigma spin up 6 3 0 0 1 3 sigma spin down 7 4 0 1 0 1 pi(lambda=+1) spin up 8 4 0 1 1 1 pi(lambda=+1) spin down 9 4 0 -1 0 1 pi(lambda=-1) spin up 10 4 0 -1 1 1 pi(lambda=-1) spin down
Definition at line 372 of file congen_projection.f90.
References congen_data::nftw.
Referenced by projec().
|
private |
Reformats (packs) the CSF expression into the style used throughout the rest of Alchemy, that is as a set of replacements from the reference determinant. Adds this to the end of the array ntdo() from location "ndo".
On entry to this routine we have the CSF defined for us as follows (from the projection step):
The above information is complemented by the analysis in the calling routine which classifies spin orbitals in this CSF wrt the reference determinant:
So, given all of the above information, the determinants in "mdo" are processed and each expressed in the format
The output is placed into array "ndto". The length available in "ndto" is passed into the routing in "len_ndto" and this is monitored to be sure we do not overflow it.
During the process, it may be necessary to multiply "cdo" by -1 as we order spin orbitals. cdo() holds the coefficients associated with each determinant. These were constructed earlier in the spin projection - note we only receive the relevant "piece" of the cdo(0 array in the arglist, the bit for this CSF, not all of it for all CSFs, as is the case with "ndto()". "nftw" is the logical unit for the printer
Definition at line 655 of file congen_projection.f90.
Referenced by wfgntr().
|
private |
Definition at line 809 of file congen_projection.f90.
References congen_data::nftw.
Referenced by projec().
|
private |
Fills the array mop
with subset of specification of every determinant for the current CSF. Only spin-orbitals that are in open shells are used, as the rest does not participate in spin composition done later in prjct. The written spin-orbitals are also immediately sorted in non-descending order. The even or odd number of swaps needed for sorting is returned (as .false. and .true., respectively) via the logical array flip
, so that the sign of determinant coefficients within this CSF can be adjusted to the used order.
OUTPUT IDOP NO OF SO IN DR BUT NOT IN DC MDOP(NELT) SO IN DR BUT NOT IN DC IDCP NO OF SO IN DC BUT NOT IN DR MDCP(NELT) SO IN DC BUT NOT IN DR IELTP NO OF SO IN OPEN SHELLS FOR A DTR MOP(MOPMX) SO
/OWF/ IDOP,IDCP,IELTPwhich was used to pass three integer values back to the caling routine. these have now been placed into the argument list; correspondingly the routine WFGNTR has been modified. It is the only routine whihc calls this one. the purpose of the variables is as follows:
IDOP holds the number of entries in MDOP IDCP holds the number of entries in MDCP IELTP is the number of electrons in open shell
Definition at line 949 of file congen_projection.f90.
References qsort().
Referenced by wfgntr().
|
private |
This routine applies the Lowdin projection operator. More details can be found in the literature at for example:
Nelson F Beebe and Sten Lucil, J Phys B: At Mol Phys, Vol. 8, Issue 14, 1975, p2320
This routine is called when a CSF is found to have two, or more electrons in open shells. Each pair of spin-orbitals in each determinant is examined and potentially used to create a new determinant. Thus the output expression for the CSF
nodo, cdo(), ndo()
may be much larger than the input
nodi, cdi(), ndo()
Note reuse of ndo() here. cdi() is used an an extendable buffer too and must have more than "nodi" elements.
The generated list of determinants is examined for any with very small coefficients (thres) and these are removed. Thus the the list may grow and shrink in this routine.
Of course the nature of the projection process is controlled by the quantum numbers input.
The routine terminates with an error message if any error conditions are found.
Definition at line 1310 of file congen_projection.f90.
References cntrct(), qsort(), rfltn(), snrm2(), and stmrg().
Referenced by wfgntr().
subroutine, public congen_projection::projec | ( | character(len=80) | sname, |
integer | megul, | ||
integer | symtyp, | ||
integer | mgvn, | ||
real(kind=wp) | s, | ||
real(kind=wp) | sz, | ||
real(kind=wp) | r, | ||
real(kind=wp) | pin, | ||
integer | nocsf, | ||
integer | byproj, | ||
integer | idiag, | ||
integer, dimension(6) | npflg, | ||
real(kind=wp) | thres, | ||
integer | nelt, | ||
integer | nsym, | ||
integer, dimension(nsym) | nob, | ||
integer, dimension(nelt) | ndtrf, | ||
integer | nftw, | ||
integer | iposit, | ||
integer, dimension(nsym) | nob0, | ||
integer, dimension(*) | nob1, | ||
integer, dimension(nsym) | nob01, | ||
integer, intent(in) | iscat, | ||
integer, intent(inout) | ntgsym, | ||
integer, dimension(ntgsym) | notgt, | ||
integer, dimension(ntgsym) | nctgt, | ||
integer, dimension(ntgsym) | mcont, | ||
integer, dimension(ntgsym) | gucont, | ||
integer, dimension(ntgsym) | mrkorb, | ||
integer, dimension(ntgsym) | mdegen, | ||
integer, intent(in) | mflag, | ||
integer, dimension(nsym) | nobe, | ||
integer, dimension(nsym) | nobp, | ||
integer, dimension(nsym) | nobv, | ||
integer | maxtgsym | ||
) |
The subroutine projec
controls the projection of the wavefunctions and writes out the final wavefunctions plus header information for future use.
Definition at line 1646 of file congen_projection.f90.
References dophz(), dophz0(), iphase(), mkorbs(), pmkorbs(), ptpwf(), rdwf(), rdwf_getsize(), wfgntr(), wrnfto(), and wrwf().
Referenced by congen_driver::csfgen().
|
private |
This provides a complete text dump of all CSFs in terms of their packed determinants.
nftw | File unit for text output of the program. |
nocsf | Number of wave functions (CSFs). |
nelt | Number of electrons (and size of ndtrf ). |
ndtrf | Reference determinant (spinorbitals per electron). |
nodi | Array with number of determinants per CSF. |
indi | Array with offsets in ndi per CSF. |
icdi | Array with offsets in cdi per CSF. |
ndi | Packed determinants. |
cdi | Deteminant coefficients. |
Definition at line 2614 of file congen_projection.f90.
Referenced by projec().
|
private |
Sort array in-place in non-descending order. Currently implemented as a stable insertion sort. This algorithm is advantageous for short and almost sorted arrays, which is the case for augmented open-shell-only subsets of determinants in CONGEN, typically resulting in asymptotic complexity of O(n).
n | Length of the array to sort. |
a | Integer array to sort. |
Definition at line 2657 of file congen_projection.f90.
Referenced by popnwf(), prjct(), and rfltn().
|
private |
Routine congen_distribution::wfn stores the CSF data in buffers of fixed size.
When the buffers are full, they are emptied to disk and reused. This is why we can have several sets to read here.
Definition at line 2699 of file congen_projection.f90.
Referenced by projec().
|
private |
Reads the information giving the size of the data used for the wavefunction, for example the number of determinants, but does not read the actual data, such as the determinants.
This is really just a stripped down version of the routine rdwf which reads the full data.
iunit | The logical unit on which the wavefucntion data is located. |
num_csfs | On return, number of CSFs stored in the file. |
num_dets | On return, number of determinants summed over all CSFs in the file. |
len_dets | On return, total summed length of all arrays of packed determinants in the file. |
Definition at line 2748 of file congen_projection.f90.
Referenced by projec().
|
private |
|
private |
This is actually an in-house version of the BLAS level 1 routine of the same name, included here to avoid dependency on BLAS library.
n | Length of the vector |
array | Vector of real numbers |
istep | Stride |
Definition at line 2844 of file congen_projection.f90.
Referenced by prjct(), and wfgntr().
|
private |
Given an existing list of determinants in cdo()/ndo() and a new single determinant cdi/ndi(), the new determinant is merged into the list and the list extended if necessary.
This operation is used during spin-projection of a CSF.
nelt | Number of electrons in each det. |
maxcdo | Dimension of cdo. |
maxndo | Dimension of ndo. |
ndo | List of determinants each with "nelt" spin-orbitals. |
cdo | Coefficient for each det in ndo. |
nodo | On input is the number of determinants in cdo/ndo. Will be updated for output if the data in cdi/ndi new entry. |
ndi | A single det of "nelt" spin orbs which has to be merged into ndo. |
cdi | Single coefficient going with the single determinant defined in ndi. |
bst | Binary search tree used for fast localization of determinants. |
Definition at line 2888 of file congen_projection.f90.
Referenced by prjct(), and rfltn().
|
private |
Takes the wavefunction generated by CSFGEN and transforms it to be fully in accord with the spin quantum numbers of the system.
Definition at line 3003 of file congen_projection.f90.
References pkwf(), popnwf(), prjct(), and snrm2().
Referenced by projec().
|
private |
Definition at line 3615 of file congen_projection.f90.
Referenced by projec().
|
private |
This subroutine writes all determinants in the current CSF to a binary file that has the following format:
[WP] n1, [INTEGER] nodo(1:n1) [WP] n2, [WP] cdo(1:n2) [WP] n3, [INTEGER] ndo(1:n3)
On entry, the file will be rewindws to its beginning. On return, the output file is rewinded to its beginning, again.
nft | An open binary file for output. |
n1 | Length of nodo (number of determinants). |
nodo | Determinant sizes. |
n2 | Length of cdo (number of determinants). |
cdo | Determinant factors within the CSF. |
n3 | Length of ndo (number of integers defining the determinants). |
ndo | Packed determinants. |
Definition at line 3739 of file congen_projection.f90.
Referenced by projec().