fock.F90 7.25 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
!!-------------------------------------------------------
!!---- Relaxed Selected Excitation (RelaxSE)
!!-------------------------------------------------------
!!---- This file is part of RelaxSE
!!---- 
!!---- The RelaxSE project is distributed under LGPL. In agreement with the
!!---- Intergovernmental Convention of the ILL, this software cannot be used
!!---- in military applications.
!!---- 
!!---- Copyright (C) 2016-2021  Institut Laue-Langevin (ILL), Grenoble, FRANCE
!!----                          Institut Neel - CNRS-UPR2940 (CNRS), Grenoble, FRANCE
!!---- 
!!---- Authors: Elisa REBOLINI (ILL)             rebolini@ill.fr
!!----          Marie-Bernadette LEPETIT (CNRS)  Marie-Bernadette.Lepetit@neel.cnrs.fr
!!---- 
!!---- RelaxSE is free software; you can redistribute it and/or
!!---- modify it under the terms of the GNU Lesser General Public
!!---- License as published by the Free Software Foundation; either
!!---- version 3.0 of the License, or (at your option) any later version.
!!---- 
!!---- RelaxSE is distributed in the hope that it will be useful,
!!---- but WITHOUT ANY WARRANTY; without even the implied warranty of
!!---- MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
!!---- Lesser General Public License for more details.
!!---- 
!!---- You should have received a copy of the GNU Lesser General Public
!!---- License along with this library; if not, see <http://www.gnu.org/licenses/>.
!!---- 


module fockmatrix

  use info
  use utils_twoint
  use detact

  implicit none

contains
  
  !$====================================================================
  !> @brief Add the 2e- integrals to the bare h_pq to build f_pq
  !> @author ER
  !> @date April 2018
  !$==================================================================== 
  subroutine build_fock(fock, hcoeur, o_info, int_info, prog_info)
    
    Real(KIND=kd_dble), dimension(:,:), allocatable, intent(inout) :: fock
    Real(KIND=kd_dble), dimension(:,:), allocatable, intent(in) :: hcoeur
    type(o_infotype), intent(in)         :: o_info
    type(int_infotype), intent(in)       :: int_info
    type(prog_infotype), intent(in)      :: prog_info

    type(intblock) :: a, ax
    integer :: i,j,p
    integer :: no,na,nv,ng

    ng = o_info%ngel
    no = o_info%nocc + o_info%nligo
    na = o_info%nact
    nv = o_info%nvirt + o_info%nligv

    fock(:,:) = hcoeur(:,:)
    
    !OO (oo|oo)
    call get_twoint(a,'oooo',o_info, int_info, prog_info%id_cpu)
    do i=ng+1,ng+no
       do j=i,ng+no
          do p = ng+1,ng+no
             fock(i,j) = fock(i,j) + 2*ijkl(a,i,j,p,p) - ijkl(a,i,p,j,p)
          enddo
          fock(j,i) = fock(i,j)
       enddo
    enddo 
    call intblock_free(a)
    
    !OA (ao|oo)
    call get_twoint(a,'aooo',o_info, int_info, prog_info%id_cpu)
    do i=ng+no+1,ng+no+na
       do j=ng+1,ng+no
          do p = ng+1,ng+no
             fock(i,j) = fock(i,j) + 2*ijkl(a,i,j,p,p) - ijkl(a,i,p,j,p)
          enddo
          fock(j,i) = fock(i,j)
       enddo
    enddo
    call intblock_free(a)

    !OV (vo|oo)
    call get_twoint(a,'vooo',o_info, int_info, prog_info%id_cpu)
    do i=ng+no+na+1,ng+no+na+nv
       do j=ng+1,ng+no
          do p = ng+1,ng+no
             fock(i,j) = fock(i,j) + 2*ijkl(a,i,j,p,p) - ijkl(a,i,p,j,p)
          enddo
          fock(j,i) = fock(i,j)
       enddo
    enddo 
    call intblock_free(a)
    
     
    !AA (aa|oo) (ao|ao)
    call get_twoint(a,'aaoo',o_info, int_info, prog_info%id_cpu)
    call get_twoint(ax,'aoao',o_info, int_info, prog_info%id_cpu)
    do i=ng+no+1,ng+no+na
       do j=i,ng+no+na
          do p = ng+1,ng+no
             fock(i,j) = fock(i,j) + 2*ijkl(a,i,j,p,p) - ijkl(ax,i,p,j,p)
          enddo
          fock(j,i) = fock(i,j)
       enddo
    enddo 
    call intblock_free(a)
    call intblock_free(ax)
    
    !AV (va|oo) (vo|ao)
    call get_twoint(a,'vaoo',o_info, int_info, prog_info%id_cpu)
    call get_twoint(ax,'voao',o_info, int_info, prog_info%id_cpu)
    do i=ng+no+na+1,ng+no+na+nv
       do j=ng+no+1,ng+no+na
          do p = ng+1,ng+no
             fock(i,j) = fock(i,j) + 2*ijkl(a,i,j,p,p) - ijkl(ax,i,p,j,p)
          enddo
          fock(j,i) = fock(i,j)
       enddo
    enddo 
    call intblock_free(a)
    call intblock_free(ax)
    
    !VV (vv|oo) (vo|vo)
    call get_twoint(a,'vvoo',o_info, int_info, prog_info%id_cpu)
    call get_twoint(ax,'vovo',o_info, int_info, prog_info%id_cpu)
    do i=ng+no+na+1,ng+no+na+nv
       do j=i,ng+no+na+nv
          do p = ng+1,ng+no
             fock(i,j) = fock(i,j) + 2*ijkl(a,i,j,p,p) - ijkl(ax,i,p,j,p)
          enddo
          fock(j,i) = fock(i,j)
       enddo
    enddo 
    call intblock_free(a)
    call intblock_free(ax)

    if ((prog_info%id_cpu .eq. 0) .and. (prog_info%iprint .gt. 1)) then
145
146
       do j=i,ng+no+na+nv
          do i=ng+1,ng+no+na+nv
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
             write(f_fock,*) i,j,fock(i,j)
          enddo
       enddo
    endif
  end subroutine build_fock

  !$====================================================================
  !> @brief Diagonalise the Fock matrix
  !> @author Elisa Rebolini
  !> @date Apr 2018
  !$====================================================================
  subroutine diag_fock(fock, o_info)
    real(kd_dble), dimension(:,:), allocatable :: fock
    type(o_infotype), intent(in) :: o_info

    logical :: debug = .false.
    integer :: n, lwork, info,i
    real(kd_dble), allocatable :: a(:,:)
    real(kd_dble), allocatable :: wr(:), work(:)
    
    n=o_info%ntot
    lwork = 3*n
    allocate(a(n,n))
    allocate(wr(n))
    allocate(work(lwork))
    a(:,:) = fock(:,:)
     
    call dsyev ('N', 'U', n, a, n, wr, work, lwork, info)
    if (debug) then
       write(f_output,*) 'Eigenvalues of the Fock Matrix'
       do i = 1,n 
          write(f_output,*) i, wr(i)
       enddo
    endif
    
    deallocate(a)
    deallocate(wr)
    deallocate(work)
    
  end subroutine diag_fock

  !$====================================================================
  !> @brief Compute E0
  !> @author Elisa Rebolini
  !> @date Apr 2018
  !$====================================================================
  subroutine get_e0(ener_info, hcoeur, o_info, int_info, prog_info)
    
    type(ener_infotype), intent(inout) :: ener_info
    real(kd_dble), dimension(:,:), allocatable, intent(in) :: hcoeur
    type(o_infotype), intent(in) :: o_info
    type(int_infotype), intent(in) :: int_info
    type(prog_infotype), intent(in) :: prog_info

    integer :: no, i, j, ng
    type(intblock) :: a
    real(kd_dble) :: e1,e2

    e1 = 0.d0
    e2 = 0.d0
    no = o_info%nocc + o_info%nligo
    ng = o_info%ngel
    call get_twoint(a,'oooo',o_info, int_info, prog_info%id_cpu)

    do i = ng+1, ng+no
       e1 = e1 + 2*hcoeur(i,i)
    enddo 
    
    do i = ng+1,ng+no
       do j = ng+1,ng+no        
             e2 = e2 + 2*ijkl(a,i,i,j,j) - ijkl(a,i,j,i,j)
       enddo
    enddo
    
    ener_info%Ecoeur = e1+e2

    if (prog_info%id_cpu.eq.0) then
       write(f_output,*) ' One-electron energy = ',e1
       write(f_output,*) ' Two-electron energy = ',e2
       write(f_output,*) ' E0                  = ',ener_info%Ecoeur
       write(f_output,*) ' Pot. Nuc.           = ',ener_info%potnuc
       write(f_output,*) ' Total E             = ',&
            ener_info%potnuc + ener_info%Ecoeur
    endif
    call intblock_free(a)
    
  end subroutine get_e0

  
end module fockmatrix