GBTOlib: library for evaluation of molecular integrals in mixed Gaussian / B-spline basis 111
Loading...
Searching...
No Matches
symmetry_gbl Module Reference

Data Types

module  geometry_obj
 Input data for the symmetry_objinit routine. More...
module  symmetry_obj
 This structure holds all necessary symmetry information on the symmetry of the molecule that can be used in the integral calculation. More...

Functions/Subroutines

integer function check_geometry_obj (this)
integer function, public determine_pg (no_sym_op, sym_op, print_to_stdout)
 Find point group by symmetry operations.
integer function, public Xlm_mgvn (n, s, l, m)
 Get (zero-based) irreducible representation of a real spherical harmonic.
logical function, public is_centre_symmetric (centre, no_sym_op, sym_op)
 Check if centre is symmetric w.r.t the symmetry operations.

Function/Subroutine Documentation

◆ check_geometry_obj()

integer function symmetry_gbl::check_geometry_obj ( class(geometry_obj) this)
Here is the call graph for this function:

◆ determine_pg()

integer function, public symmetry_gbl::determine_pg ( integer no_sym_op,
character(len=sym_op_nam_len), dimension(3) sym_op,
logical, optional print_to_stdout )

Find point group by symmetry operations.

Authors
Zdenek Masin
Date
2016

Returns the group identifier as defined in "const.f90", based on the set of symmetry operations passed as arguments (each one of 'X', 'Y', 'Z', 'XY', 'YZ', 'XZ', 'XYZ', indicating simple and combined plane reflections that produce no extra sign).

Parameters
no_sym_opNumber of symmetry operations set in sym_op.
sym_opA triplet of strings, of which only no_sym_op need to be set.
print_to_stdoutPrint the group name to stdout.
Here is the call graph for this function:

◆ is_centre_symmetric()

logical function, public symmetry_gbl::is_centre_symmetric ( real(cfp), dimension(3), intent(in) centre,
integer, intent(in) no_sym_op,
character(len=sym_op_nam_len), dimension(:), intent(in) sym_op )

Check if centre is symmetric w.r.t the symmetry operations.

Authors
J Benda
Date
2025

Given that all symmetry operations in UKRmol+ are composed of axis inversions, the check needs to verify simply that whenever a specific axis inversion is included in any composed symmetry operation, the corresponding component of centre position is zero.

◆ Xlm_mgvn()

integer function, public symmetry_gbl::Xlm_mgvn ( integer, intent(in) n,
character(len=sym_op_nam_len), dimension(:), intent(in) s,
integer, intent(in) l,
integer, intent(in) m )

Get (zero-based) irreducible representation of a real spherical harmonic.

Authors
J Benda
Date
2025

We use the symmetry properties of the real spherical harmonics. Assuming that X, Y and Z are mirror operations by coordinate planes orthogonal to the given axes, i.e.,

  • X: phi -> pi - phi
  • Y: phi -> -phi
  • Z: theta -> pi - theta

the symmetry properties of the real spherical harmonics are:

  • X X_{l,-|m|} = -(-1)^m X_{l,-|m|}
  • Y X_{l,-|m|} = -X_{l,+|m|}
  • Z X_{l,-|m|} = (-1)^(l-m) X_{l,-|m|}
  • X X_{l,0} = X_{l,0}
  • Y X_{l,0} = X_{l,0}
  • Z X_{l,0} = (-1)^(l-m) X_{l,0}
  • X X_{l,+|m|} = (-1)^m X_{l,+|m|}
  • Y X_{l,+|m|} = X_{l,+|m|}
  • Z X_{l,+|m|} = (-1)^(l-m) * X_{l,+|m|}

Of these, only the specified symmetry elements are used. The relevant characters of the supported groups are the following (consistent with Molpro convention):

  • C1: A
  • Ci: Ag (XYZ+), Au (XYZ-)
  • Cs: A' (X+), A'' (X-)
  • C2: A (XY+), B (XY-)
  • D2: A (XZ+,YZ+), B3 (XZ-,YZ+), B2 (XZ+,YZ-), B1 (XZ-,YZ-)
  • C2h: Ag (Z+,XY+), Au (Z-,XY+), Bu (Z+,XY-), Bg (Z-,XY-)
  • C2v: A (X+Y+), B1 (X-Y+), B2 (X+Y-), A2 (X-Y-)
  • D2h: Ag (X+Y+Z+), B3u (X-Y+Z+), B2u (X+Y-Z+), B1g (X-Y-Z+), B1u (X+Y+Z-), B2g (X-Y+Z-), B3g (X+Y-Z-), Au (X-Y-Z-)

To maintain a predictable pattern in these tables, the symmetry operations have to be sorted by length (i.e. 'Z' first, then 'XY') and by alphabet ('X' first, then 'Y').