MPI-SCATCI  2.0
An MPI version of SCATCI
postprocessing_module Module Reference

Further processing of the diagonalization results. More...

Data Types

type  outerinterface
 SWINTERF-inspired post-processing object. More...
 

Functions/Subroutines

subroutine, public postprocess (SCATCI_input, solutions)
 Full post-processing pass. More...
 
subroutine redistribute_solutions (SCATCI_input, solutions)
 Redistribute solutions from groups to everyone. More...
 
subroutine cdenprop_properties (SCATCI_input, solutions, properties_all)
 Evaluate multipoles in CDENPROP. More...
 
subroutine outer_interface (SCATCI_input, solutions, inner_properties)
 Extract data needed by outer-region codes. More...
 
subroutine init_outer_interface (this, input)
 Allocate memory for channel information. More...
 
subroutine deinit_outer_interface (this)
 Release memory held by this object. More...
 
subroutine setup_amplitudes (this, opts)
 Initialize raw boundary amplitudes (in integral library) More...
 
subroutine extract_data (this, opts, prop, solution)
 Get interface data. More...
 
subroutine write_data (this, i, opts, prop, solution)
 Write interface data. More...
 
subroutine get_channel_info (this, opts, prop)
 Assemble the list of outer channels. More...
 
subroutine get_boundary_data (this, opts, prop, solution)
 Evaluate boundary amplitudes for propagation and RMT. More...
 
subroutine write_channel_info (this, nchset, opts, cform, prop)
 Write the channel list to disk. More...
 
subroutine write_boundary_data (this, nrmset, opts, rform, prop, solution)
 Write the R-matrix amplitudes to disk. More...
 
subroutine get_channel_couplings (this, ismax, ntarg, target_properties, alpha0, alpha2, use_pol)
 Evaluate long-range channel coupling coefficients. More...
 
subroutine write_rmt_data (input, inner_properties, target_properties, solutions, intf)
 Compose the RMT molecular input data file. More...
 
subroutine generate_couplings (maxl, n_rg, rg, lm_rg)
 Evaluate angular couplings for outer region of RMT. More...
 

Variables

integer, parameter rmt_int = shortint
 
integer, parameter rmt_real = wp
 

Detailed Description

Further processing of the diagonalization results.

Authors
J Benda
Date
2019

Contains subroutines that can be used to extract further interesting data from the eigenvectors. In particular, this means extraction of boundary amplitudes for outer region codes (see outer_interface) and calculation of multipole moments, mainly for use in the time-dependent R-matrix code (see cdenprop_properties).

This module currently supports only the UKRmol+ integral format.

Function/Subroutine Documentation

◆ cdenprop_properties()

subroutine postprocessing_module::cdenprop_properties ( type(optionsset), intent(in)  SCATCI_input,
type(solutionhandler), dimension(:), intent(inout), allocatable  solutions,
type(molecular_properties_data), intent(inout)  properties_all 
)

Evaluate multipoles in CDENPROP.

Authors
J Benda
Date
2019

Calculate properties (multipole moments) of the solutions using CDENPROP library. Only those symmetries that have igh = 1, vecstore = 3 set will be used. This subroutine also assumes that setups for all symmetries are compatible, in particular:

  • use the same molecular integral format
  • use the same CI targets (in case of scattering diagonalization)
Parameters
[in]SCATCI_inputInput SCATCI namelist data for all symmetries.
[in]solutionsEigenvectors and related stuff calculated by the diagonalizer.
[out]properties_allTable of properties calculated by CDENPROP.

Definition at line 225 of file Postprocessing_module.F90.

Here is the caller graph for this function:

◆ deinit_outer_interface()

subroutine postprocessing_module::deinit_outer_interface ( class(outerinterface), intent(inout)  this)

Release memory held by this object.

Authors
J Benda
Date
2019
Parameters
[in]thisInterface object to finalize.

Definition at line 459 of file Postprocessing_module.F90.

◆ extract_data()

subroutine postprocessing_module::extract_data ( class(outerinterface), intent(inout)  this,
type(options), intent(in)  opts,
type(molecular_properties_data), intent(in)  prop,
type(solutionhandler), intent(in)  solution 
)

Get interface data.

Authors
J Benda
Date
2019

This mirros the functionality of SWITERF (UKRmol-out), particularly:

  • reads the target property file (fort.24)
  • assembles the list of outer region channels
  • evaluates amplitudes of the eigenstates at R-matrix boundary and, if requested, in its vicinity for use in RMT

All data are then stored internally for writing.

Parameters
[in]thisInterface object to read.
[in]optsSCATCI input namelist data.
[in]propTarget properties.
[in]solutionEigenvectors and related stuff calculated by the diagonalizer.

Definition at line 534 of file Postprocessing_module.F90.

◆ generate_couplings()

subroutine postprocessing_module::generate_couplings ( integer(rmt_int), intent(in)  maxl,
integer(rmt_int), intent(out)  n_rg,
real(rmt_real), dimension(:), intent(inout), allocatable  rg,
integer(rmt_int), dimension(:,:), intent(inout), allocatable  lm_rg 
)

Evaluate angular couplings for outer region of RMT.

Authors
Z Masin
Date
2017

Precompute the necessary Gaunt coefficients.

Parameters
[in]maxlLimit on angular quantum number.
[out]n_rgNumber of calculated Gaunt coefficients.
[out]rgValues of the non-zero Gaunt coefficients.
[out]lm_rgInteger sextets l1,m1,l2,m2,l3,m3 for each of the non-zero Gaunt coefficients.

Definition at line 1694 of file Postprocessing_module.F90.

Here is the caller graph for this function:

◆ get_boundary_data()

subroutine postprocessing_module::get_boundary_data ( class(outerinterface), intent(inout)  this,
type(options), intent(in)  opts,
type(molecular_properties_data), intent(in)  prop,
type(solutionhandler), intent(in)  solution 
)

Evaluate boundary amplitudes for propagation and RMT.

Authors
J Benda
Date
2019

Calculates the contribution to each channel amplitude per eigenstate. The master process retrieves the raw boundary amplitudes from the integral library, then distributes them in the BLACS context of the solution vector, where the multiplication takes place.

The boundary amplitudes evaluated at the R_matrix sphere (but not those evaluated inside) are then collected to master process, who will need to write them to the record-based swinterf output file.

Parameters
[in]thisInterface object to update.
[in]optsSCATCI namelist information for this irreducible representation.
[in]propResults from denprop.
[in]solutionEigenvectors and related stuff calculated by the diagonalizer.

Definition at line 734 of file Postprocessing_module.F90.

◆ get_channel_couplings()

subroutine postprocessing_module::get_channel_couplings ( class(outerinterface), intent(inout)  this,
integer, intent(in)  ismax,
integer, intent(in)  ntarg,
type(molecular_properties_data), intent(in)  target_properties,
real(wp), intent(in)  alpha0,
real(wp), intent(in)  alpha2,
logical, intent(in)  use_pol 
)

Evaluate long-range channel coupling coefficients.

Authors
Z Masin
Date
2014

This is a completely rewritten version of the original SWINTERF routine SWASYMC. The coupling potentials for the outer region calculation are defined as:

\[ V_{ij}(r)=\sum_{\lambda=0}^{\infty}\frac{1}{r^{\lambda+1}}\times \nonumber \\ \times\underbrace{\sum_{m=-\lambda}^{\lambda}\langle{\cal Y}_{l_{i},m_{i}}\vert {\cal {Y}}_{\lambda ,m}\vert{\cal Y}_{l_{j},m_{j}}\rangle \sqrt{\frac{4\pi}{2\lambda+1}}\underbrace{\sqrt{\frac{4\pi}{2\lambda+1}}\left( T_{ij}^{\lambda m} - \langle \Phi_{i}\vert\Phi_{j}\rangle \sum_{k=1}^{Nuclei}Z_{k}{\cal Y}_{\lambda,m}(\hat{\mathbf{R}}_{k})R_{k}^{\lambda}\right)}_{Q_{ij}^{\lambda m}}}_{a_{ij\lambda}} . \]

This routine calculates the coefficients \( a_{ij\lambda} \) using the information for all channels (i,j) and the target permanent and transition multipole moments \( Q_{ij}^{\lambda m} \). The main difference from SWASYMC is that here the coupling coefficients for the real spherical harmonics are calculated independently of the molecular orientation, symmetry and the lambda value. Tests were performed for pyrazine ( \( D_{2h} \)), uracil ( \( C_s \)) and water ( \( C_{2v} \)). The inclusion of the additional polarization potential is also possible, but note that only the spherical part is taken into account at the moment (see the comment for the variable alpha2).

Note
06/02/2019 - Jakub Benda: Copied over from UKRmol-out.
Parameters
[in,out]thisInterface object to update.
[in]ismaxMaximum lambda value for the multipole potential contribution.
[in]ntargTotal number of target states.
[in]target_propertiesDenprop data.
[in]alpha0Spherical part of the ground state target polarizability.
[in]alpha2Non-spherical part of the ground state target polarizability. Note that this has not been implemented in this subroutine yet since it is not clear to me what is the convention defining this value in the old code.
[in]use_polIf .true. then the coefficients for the polarization potential will be constructed.

Definition at line 1140 of file Postprocessing_module.F90.

◆ get_channel_info()

subroutine postprocessing_module::get_channel_info ( class(outerinterface), intent(inout)  this,
type(options), intent(in)  opts,
type(molecular_properties_data), intent(in)  prop 
)

Assemble the list of outer channels.

Authors
J Benda
Date
2019

Adapted from SWCHANL (UKRmol-out). Reads the target information from fort.24 and, with the knowledge of the (N+1)-electron setup, determines the individual channels.

Parameters
[in,out]thisInterface object to update.
[in]optsSCATCI namelist information for this irreducible representation.
[in]propTarget properties.

Definition at line 616 of file Postprocessing_module.F90.

◆ init_outer_interface()

subroutine postprocessing_module::init_outer_interface ( class(outerinterface), intent(inout)  this,
type(optionsset), intent(in)  input 
)

Allocate memory for channel information.

Authors
J Benda
Date
2019
Parameters
[in]thisInterface object to initialize.
[in]inputInput SCATCI namelist data for all symmetries.

Definition at line 410 of file Postprocessing_module.F90.

◆ outer_interface()

subroutine postprocessing_module::outer_interface ( type(optionsset), intent(in)  SCATCI_input,
type(solutionhandler), dimension(:), intent(in), allocatable  solutions,
type(molecular_properties_data), intent(in)  inner_properties 
)

Extract data needed by outer-region codes.

Authors
J Benda
Date
2019

This subroutine does the following:

  • write channel information to unit LUCHAN
  • write boundary amplitudes to unit LURMT
  • write RMT data to file molecular_data

Essentially, it should cover the work of SWINTERF and RMT_INTERFACE. At the moment, only the symmetries marked with vecstore = 3 are processed in the outer interface.

Parameters
[in]SCATCI_inputInput SCATCI namelist data for all symmetries.
[in]solutionsEigenvectors and related stuff calculated by the diagonalizer.
[out]inner_propertiesTable of (N+1) properties calculated by CDENPROP.

Definition at line 357 of file Postprocessing_module.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ postprocess()

subroutine, public postprocessing_module::postprocess ( type(optionsset), intent(in)  SCATCI_input,
type(solutionhandler), dimension(:), intent(inout), allocatable  solutions 
)

Full post-processing pass.

Authors
J Benda
Date
2019

Call the post-processing sequence, starting with evaluation of inner-region dipole moments, followed by extraction of information useful for the outer region.

Parameters
[in]SCATCI_inputInput SCATCI namelist data for all symmetries.
[in]solutionsEigenvectors and related stuff calculated by the diagonalizer.

Definition at line 103 of file Postprocessing_module.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ redistribute_solutions()

subroutine postprocessing_module::redistribute_solutions ( type(optionsset), intent(in)  SCATCI_input,
type(solutionhandler), dimension(:), intent(inout), allocatable  solutions 
)

Redistribute solutions from groups to everyone.

Authors
J Benda
Date
2019

This subroutine needs to be called before using solutions in other subroutines of this module. The diagonalizations may have run in MPI groups and so the eigenvectors are not evenly distributed among all processes, which would cause comlications in communication. This subroutine will redistribute the eigenvectors from groups to the MPI world group.

Parameters
[in]SCATCI_inputInput SCATCI namelist data for all symmetries.
[in]solutionsEigenvectors and related stuff calculated by the diagonalizer.

Definition at line 132 of file Postprocessing_module.F90.

Here is the caller graph for this function:

◆ setup_amplitudes()

subroutine postprocessing_module::setup_amplitudes ( class(outerinterface), intent(inout)  this,
type(options), intent(in)  opts 
)
private

Initialize raw boundary amplitudes (in integral library)

Authors
J Benda
Date
2019
Parameters
[in]thisInterface object to read.
[in]optsSCATCI input namelist data.

Definition at line 502 of file Postprocessing_module.F90.

◆ write_boundary_data()

subroutine postprocessing_module::write_boundary_data ( class(outerinterface), intent(in)  this,
integer, intent(in)  nrmset,
type(options), intent(in)  opts,
character(len=1), intent(in)  rform,
type(molecular_properties_data), intent(in)  prop,
type(solutionhandler), intent(in)  solution 
)

Write the R-matrix amplitudes to disk.

Authors
J Benda
Date
2019

Writes the information into the selected set of the R-matrix output file (fort.21). Adapted from WRITRM (UKRmol-out).

Parameters
[in]thisInterface object to update.
[in]nrmsetOutput data set in the file unit.
[in]optsSCATCI namelist information for this irreducible representation.
[in]rform'F'/'U' for formatted or unformatted output.
[in]propResults from denprop.
[in]solutionEigenvectors and related stuff calculated by the diagonalizer.
Todo:
At the moment, this subroutine is called (and boundary amplitudes written) by master processs only, which holds all boundary amplitudes. In cases with many channels and many inner region eigenstates, the size of the boundary amplitude matrices can, hypothetically, become too large for a single process to manage (either because of memory or I/O bandwidth limitations). It will then be necessary to keep the amplitudes distributed (in get_boundary_data) and write them using MPI I/O similarly as in write_rmt_data. Note that in this case we are writing a SEQUENTIAL rather than a STREAM file, so the binary data written need to be split into chunks of size at most 2 GiB, each of them immediately preceded by a "bookmark" 32-bit integer with the length of the chunk in bytes, and followed by a similar bookmark with the same value, but with positive sign when that chunk is the last one in the record, or negative sign when that chunk is yet followed by another chunk within the same record. While not elaborated in the Fortran standard, this seems to be the universal storage scheme used in SEQUENTIAL UNFORMATTED files by GNU, Intel and Cray compilers.

Definition at line 1007 of file Postprocessing_module.F90.

◆ write_channel_info()

subroutine postprocessing_module::write_channel_info ( class(outerinterface), intent(in)  this,
integer, intent(in)  nchset,
type(options), intent(in)  opts,
character(len=1), intent(in)  cform,
type(molecular_properties_data), intent(in)  prop 
)

Write the channel list to disk.

Authors
J Benda
Date
2019

Adapted from WRITCH (UKRmol-out). Writes the information into the selected set of the channel output file fort.10.

Parameters
[in]thisInterface object to update.
[in]nchsetOutput data set in the file unit.
[in]optsSCATCI namelist information for this irreducible representation.
[in]cform'F'/'U' for formatted or unformatted output.
[in]propResults from denprop.

Definition at line 888 of file Postprocessing_module.F90.

◆ write_data()

subroutine postprocessing_module::write_data ( class(outerinterface), intent(in)  this,
integer, intent(in)  i,
type(optionsset), intent(in)  opts,
type(molecular_properties_data), intent(in)  prop,
type(solutionhandler), intent(in)  solution 
)

Write interface data.

Authors
J Benda
Date
2019

Updates the channel and R-matrix dump files (fort.10 and fort.21, respectively) by adding/overwriting a data set.

Parameters
[in]thisInterface object to read.
[in]iOutput set number for this irreducible representation.
[in]optsSCATCI input namelist data.
[in]propTarget properties.
[in]solutionEigenvectors and related stuff calculated by the diagonalizer.

Definition at line 584 of file Postprocessing_module.F90.

◆ write_rmt_data()

subroutine postprocessing_module::write_rmt_data ( type(optionsset), intent(in)  input,
type(molecular_properties_data), intent(in)  inner_properties,
type(molecular_properties_data), intent(in)  target_properties,
type(solutionhandler), dimension(:), intent(in), allocatable  solutions,
type(outerinterface), dimension(:), intent(in), allocatable  intf 
)

Compose the RMT molecular input data file.

Authors
J Benda
Date
2019

Based on generate_data_for_rmt from the original program rmt_interface (UKRmol-out).

This subroutine writes the input data file used by molecular version of RMT. The file is always called 'molecular_data', it is a binary stream file and it contains the following data:

  • (For m = -1, 0, 1) s, s1, s2, iidip(1:s), ifdip(1:s), dipsto(1:s1,1:s2, 1:s): Number of distinct symmetry pairs, maximal number of states in bra symmetry maximal number of states in ket symmetry, list of bra symmetries, list of ket symmetries, (N + 1)-electron state transition dipole blocks (padded by zeros to s1 x s2).
  • ntarg: Number of targets (N-electron ions).
  • (For m = -1, 0, 1) crlv: Dipole transition matrix for N-electron states
  • n_rg, rg, lm_rg: Number of real Gaunt coeffcients, array of them, their angular momentum quantum numbers.
  • nelc, nz, lrang2, lamax, ntarg, inast, nchmx, nstmx, lmaxp1: Number of electrons, nuclear charge, maximal angular momentum, maximal angular momentum transfer, number of target states, number of irreducible representations, maximal number of channels in all irreducible representations, maximal number of eigenstates in all irreducible representations, maximal angular momentum.
  • rmatr, bbloch: R-matrix radius, Bloch b-coeffcient.
  • etarg, ltarg, starg: Ionic state eigen-energies, their irreducible representations (one-based), and spins.
  • nfdm, delta_r: Number of finite difference points inside the R-matrix sphere, their spacing.
  • r_points: List of the finite difference points inside the R-matrix sphere.
  • (For all target symmetries) lrgl, nspn, npty, nchan, mnp1, nconat, l2p, m2p, eig, wamp, cf, ichl, s1, s2, wamp2: Irreducible representation (one-based), spin multiplicity, parity (not used, always zero), number of continuum channels, number of (N + 1)-electron states, number of continuum channels per target (ionic) state, angular momentum per channel, angular momentum projection per channel, eigenenergies of the (N + 1)-electron states (+ core energy), contributions of the (N + 1)-electron eigenstates to every channel (boundary amplitudes), long-range multipole coupling coefficients, target state index per continuum channel, row size of wamp2, col size of wamp2, wamp evaluated at inner-region finite difference points.
Todo:
At the moment, the subroutine assumes that all inner states have the same spin.
Warning
The subroutine makes use of MPI I/O for writing the distributed arrays. In the present, straightforward implementation this limits the local portion of the distributed matrices to 2**31 elements (i.e. 16 GiB for 8-byte real arrays). If this became a problem, one would need to rewrite the subroutine to write everything in one go as a custom MPI datatype, whose building routines accept long integers. However, 16 GiB per core per distributed matrix seems quite acceptable for now.
Parameters
[in]inputInput SCATCI namelist data for all symmetries.
[in]inner_propertiesResults from cdenprop.
[in]target_propertiesResults from denprop.
[in]solutionsResults from diagonalization.
[in]intfExtracted channel and boundary data.

Definition at line 1374 of file Postprocessing_module.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ rmt_int

integer, parameter postprocessing_module::rmt_int = shortint
private

Definition at line 45 of file Postprocessing_module.F90.

◆ rmt_real

integer, parameter postprocessing_module::rmt_real = wp
private

Definition at line 46 of file Postprocessing_module.F90.