MPI-SCATCI
2.0
An MPI version of SCATCI
|
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 |
Further processing of the diagonalization results.
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.
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.
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:
[in] | SCATCI_input | Input SCATCI namelist data for all symmetries. |
[in] | solutions | Eigenvectors and related stuff calculated by the diagonalizer. |
[out] | properties_all | Table of properties calculated by CDENPROP. |
Definition at line 225 of file Postprocessing_module.F90.
subroutine postprocessing_module::deinit_outer_interface | ( | class(outerinterface), intent(inout) | this | ) |
Release memory held by this object.
[in] | this | Interface object to finalize. |
Definition at line 459 of file Postprocessing_module.F90.
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.
This mirros the functionality of SWITERF (UKRmol-out), particularly:
All data are then stored internally for writing.
[in] | this | Interface object to read. |
[in] | opts | SCATCI input namelist data. |
[in] | prop | Target properties. |
[in] | solution | Eigenvectors and related stuff calculated by the diagonalizer. |
Definition at line 534 of file Postprocessing_module.F90.
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.
Precompute the necessary Gaunt coefficients.
[in] | maxl | Limit on angular quantum number. |
[out] | n_rg | Number of calculated Gaunt coefficients. |
[out] | rg | Values of the non-zero Gaunt coefficients. |
[out] | lm_rg | Integer sextets l1,m1,l2,m2,l3,m3 for each of the non-zero Gaunt coefficients. |
Definition at line 1694 of file Postprocessing_module.F90.
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.
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.
[in] | this | Interface object to update. |
[in] | opts | SCATCI namelist information for this irreducible representation. |
[in] | prop | Results from denprop. |
[in] | solution | Eigenvectors and related stuff calculated by the diagonalizer. |
Definition at line 734 of file Postprocessing_module.F90.
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.
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
).
[in,out] | this | Interface object to update. |
[in] | ismax | Maximum lambda value for the multipole potential contribution. |
[in] | ntarg | Total number of target states. |
[in] | target_properties | Denprop data. |
[in] | alpha0 | Spherical part of the ground state target polarizability. |
[in] | alpha2 | Non-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_pol | If .true. then the coefficients for the polarization potential will be constructed. |
Definition at line 1140 of file Postprocessing_module.F90.
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.
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.
[in,out] | this | Interface object to update. |
[in] | opts | SCATCI namelist information for this irreducible representation. |
[in] | prop | Target properties. |
Definition at line 616 of file Postprocessing_module.F90.
subroutine postprocessing_module::init_outer_interface | ( | class(outerinterface), intent(inout) | this, |
type(optionsset), intent(in) | input | ||
) |
Allocate memory for channel information.
[in] | this | Interface object to initialize. |
[in] | input | Input SCATCI namelist data for all symmetries. |
Definition at line 410 of file Postprocessing_module.F90.
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.
This subroutine does the following:
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.
[in] | SCATCI_input | Input SCATCI namelist data for all symmetries. |
[in] | solutions | Eigenvectors and related stuff calculated by the diagonalizer. |
[out] | inner_properties | Table of (N+1) properties calculated by CDENPROP. |
Definition at line 357 of file Postprocessing_module.F90.
subroutine, public postprocessing_module::postprocess | ( | type(optionsset), intent(in) | SCATCI_input, |
type(solutionhandler), dimension(:), intent(inout), allocatable | solutions | ||
) |
Full post-processing pass.
Call the post-processing sequence, starting with evaluation of inner-region dipole moments, followed by extraction of information useful for the outer region.
[in] | SCATCI_input | Input SCATCI namelist data for all symmetries. |
[in] | solutions | Eigenvectors and related stuff calculated by the diagonalizer. |
Definition at line 103 of file Postprocessing_module.F90.
subroutine postprocessing_module::redistribute_solutions | ( | type(optionsset), intent(in) | SCATCI_input, |
type(solutionhandler), dimension(:), intent(inout), allocatable | solutions | ||
) |
Redistribute solutions from groups to everyone.
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.
[in] | SCATCI_input | Input SCATCI namelist data for all symmetries. |
[in] | solutions | Eigenvectors and related stuff calculated by the diagonalizer. |
Definition at line 132 of file Postprocessing_module.F90.
|
private |
Initialize raw boundary amplitudes (in integral library)
[in] | this | Interface object to read. |
[in] | opts | SCATCI input namelist data. |
Definition at line 502 of file Postprocessing_module.F90.
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.
Writes the information into the selected set of the R-matrix output file (fort.21). Adapted from WRITRM (UKRmol-out).
[in] | this | Interface object to update. |
[in] | nrmset | Output data set in the file unit. |
[in] | opts | SCATCI namelist information for this irreducible representation. |
[in] | rform | 'F'/'U' for formatted or unformatted output. |
[in] | prop | Results from denprop. |
[in] | solution | Eigenvectors and related stuff calculated by the diagonalizer. |
Definition at line 1007 of file Postprocessing_module.F90.
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.
Adapted from WRITCH (UKRmol-out). Writes the information into the selected set of the channel output file fort.10.
[in] | this | Interface object to update. |
[in] | nchset | Output data set in the file unit. |
[in] | opts | SCATCI namelist information for this irreducible representation. |
[in] | cform | 'F'/'U' for formatted or unformatted output. |
[in] | prop | Results from denprop. |
Definition at line 888 of file Postprocessing_module.F90.
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.
Updates the channel and R-matrix dump files (fort.10 and fort.21, respectively) by adding/overwriting a data set.
[in] | this | Interface object to read. |
[in] | i | Output set number for this irreducible representation. |
[in] | opts | SCATCI input namelist data. |
[in] | prop | Target properties. |
[in] | solution | Eigenvectors and related stuff calculated by the diagonalizer. |
Definition at line 584 of file Postprocessing_module.F90.
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.
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:
[in] | input | Input SCATCI namelist data for all symmetries. |
[in] | inner_properties | Results from cdenprop. |
[in] | target_properties | Results from denprop. |
[in] | solutions | Results from diagonalization. |
[in] | intf | Extracted channel and boundary data. |
Definition at line 1374 of file Postprocessing_module.F90.
|
private |
Definition at line 45 of file Postprocessing_module.F90.
|
private |
Definition at line 46 of file Postprocessing_module.F90.