MPI-SCATCI 2.0
An MPI version of SCATCI
Loading...
Searching...
No Matches
BaseCorePotential_module.f90
Go to the documentation of this file.
1! Copyright 2019
2!
3! For a comprehensive list of the developers that contributed to these codes
4! see the UK-AMOR website.
5!
6! This file is part of UKRmol-in (UKRmol+ suite).
7!
8! UKRmol-in is free software: you can redistribute it and/or modify
9! it under the terms of the GNU General Public License as published by
10! the Free Software Foundation, either version 3 of the License, or
11! (at your option) any later version.
12!
13! UKRmol-in is distributed in the hope that it will be useful,
14! but WITHOUT ANY WARRANTY; without even the implied warranty of
15! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16! GNU General Public License for more details.
17!
18! You should have received a copy of the GNU General Public License
19! along with UKRmol-in (in source/COPYING). Alternatively, you can also visit
20! <https://www.gnu.org/licenses/>.
21
22!> \brief Base type for core potentials
23!> \authors A Al-Refaie
24!> \date 2017
25!>
26!> Used by specialized core potential types, eg. MOLPROCorePotential.
27!>
28!> \note 16/01/2019 - Jakub Benda: Unifom coding style and expanded documentation.
29!>
30module basecorepotential_module
31
32 use precisn, only: wp
33
34 implicit none
35
36 type, abstract :: basecorepotential
37 real(wp), allocatable :: target_energies(:,:)
38 real(wp) :: l2_ecp_deficit
39 integer :: num_target_sym
40 integer, allocatable :: num_target_per_symmetries(:)
41 integer, allocatable :: spatial_symmetry(:)
42 integer, allocatable :: spin_symmetry(:)
43 real(wp) :: ground_energy
44 contains
45 procedure, public :: construct => construct_core
46 procedure(generic_ecp), deferred :: parse_ecp
47 procedure(generic_input), deferred :: parse_input
48 procedure(generic_energy), deferred :: modify_target_energies
49 procedure, public :: modify_ecp_diagonal
50 procedure, public :: modify_ecp_diagonal_l2
51 end type basecorepotential
52
53 abstract interface
54 !> \brief Main build routine of the hamiltonian
55 !> \authors A Al-Refaie
56 !> \date 2017
57 !> \public
58 !>
59 !> All build must be done within this routine in order to be used by MPI-SCATCI.
60 !>
61 !> \param[out] matrix_elements Resulting matrix elements from the build
62 !>
63 subroutine generic_ecp (this, filename)
64 import :: basecorepotential
65 class(basecorepotential) :: this
66 character(len=*) :: filename
67 end subroutine generic_ecp
68 end interface
69
70 abstract interface
71 !> \brief Main build routine of the hamiltonian
72 !> \authors A Al-Refaie
73 !> \date 2017
74 !> \public
75 !>
76 !> All build must be done within this routine in order to be used by MPI-SCATCI.
77 !>
78 !> \param[out] matrix_elements Resulting matrix elements from the build
79 !>
80 subroutine generic_input(this)
81 import :: basecorepotential
82 class(basecorepotential) :: this
83 end subroutine generic_input
84 end interface
85
86 abstract interface
87 subroutine generic_energy (this, target_sym, target_energy)
88 import :: basecorepotential
89 class(basecorepotential) :: this
90 integer, intent(in) :: target_sym
91 real(wp) :: target_energy(:)
92 end subroutine generic_energy
93 end interface
94
95contains
96
97 subroutine construct_core (this, num_target_sym, num_target_per_symmetries, spatial_symmetry, spin_symmetry)
98 class(basecorepotential) :: this
99 integer, intent(in) :: num_target_sym
100 integer, intent(in) :: num_target_per_symmetries(num_target_sym), &
101 spatial_symmetry(num_target_sym), &
102 spin_symmetry(num_target_sym)
103
104 this % num_target_sym = num_target_sym
105 allocate(this % num_target_per_symmetries(num_target_sym))
106
107 this % num_target_per_symmetries = num_target_per_symmetries
108 allocate(this % target_energies(maxval(this % num_target_per_symmetries), this % num_target_sym))
109
110 allocate(this % spatial_symmetry(num_target_sym))
111 allocate(this % spin_symmetry(num_target_sym))
112
113 this % spatial_symmetry = spatial_symmetry
114 this % spin_symmetry = spin_symmetry
115
116 this % target_energies = 0.0
117 this % L2_ECP_deficit = 0.0
118 end subroutine construct_core
119
120
121 subroutine modify_ecp_diagonal (this, value)
122 class(basecorepotential), intent(in) :: this
123 real(wp), intent(inout) :: value
124
125 value = value - this % ground_energy
126 value = value + this % target_energies(1,1)
127 end subroutine modify_ecp_diagonal
128
129
130 subroutine modify_ecp_diagonal_l2 (this, value)
131 class(basecorepotential),intent(in) :: this
132 real(wp), intent(inout) :: value
133
134 value = value + this % L2_ECP_deficit
135 end subroutine modify_ecp_diagonal_l2
136
137end module basecorepotential_module