Octopus
lda_u_mixer.F90
Go to the documentation of this file.
1!! Copyright (C) 2016 N. Tancogne-Dejean
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22 use debug_oct_m
23 use global_oct_m
24 use lda_u_oct_m
26 use mix_oct_m
31 use types_oct_m
32
33 implicit none
34
35 private
36
37 public :: &
46
47 type lda_u_mixer_t
48 private
49 integer :: occsize
50 logical :: realstates
51 logical :: apply = .false.
52 logical :: mixU = .false.
53
54 real(real64), allocatable :: dtmp_occ(:,:), tmpU(:,:)
55 complex(real64), allocatable :: ztmp_occ(:,:)
56
57 type(mixfield_t) :: mixfield_occ, mixfield_U
58 end type lda_u_mixer_t
59
60contains
61
62 ! ---------------------------------------------------------
63 subroutine lda_u_mixer_init_auxmixer(this, namespace, mixer, smix, st)
64 type(lda_u_t), intent(in) :: this
65 type(namespace_t), intent(in) :: namespace
66 type(lda_u_mixer_t), intent(inout) :: mixer
67 type(mix_t), intent(inout) :: smix
68 type(states_elec_t), intent(in) :: st
69
70 integer :: dim1
71
72 if (this%level == dft_u_none) return
74
75 dim1 = this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets
76 if (this%level == dft_u_acbn0) then
77 dim1 = dim1*2
78 if (this%intersite) then
79 dim1 = dim1 + 2*this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets*this%maxneighbors
80 end if
81 end if
82
83 if (states_are_real(st)) then
84 call mixfield_init(smix, mixer%mixfield_occ, dim1, 1, mix_d3(smix), type_float)
85 mixer%realstates = .true.
86 else
87 call mixfield_init(smix, mixer%mixfield_occ, dim1, 1, mix_d3(smix), type_cmplx)
88 mixer%realstates = .false.
89 end if
90 call mixfield_clear(mix_scheme(smix), mixer%mixfield_occ)
91 call mix_add_auxmixfield(namespace, smix, mixer%mixfield_occ)
92
93 if (this%level == dft_u_acbn0) then
94 call mixfield_init(smix, mixer%mixfield_U, this%norbsets, 1, mix_d3(smix), type_float)
95 call mixfield_clear(mix_scheme(smix), mixer%mixfield_U)
96 call mix_add_auxmixfield(namespace, smix, mixer%mixfield_U)
97 mixer%mixU = .true.
98 else
99 mixer%mixU = .false.
100 end if
101
103 end subroutine lda_u_mixer_init_auxmixer
104
105
106 ! ---------------------------------------------------------
107 subroutine lda_u_mixer_init(this, mixer, st)
108 type(lda_u_t), intent(in) :: this
109 type(lda_u_mixer_t), intent(inout) :: mixer
110 type(states_elec_t), intent(in) :: st
111
112 if (this%level == dft_u_none) return
113 push_sub(lda_u_mixer_init)
115 mixer%apply = .true.
116
117 mixer%occsize = this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets
118 if (this%level == dft_u_acbn0) then
119 mixer%occsize = mixer%occsize*2
120 if (this%intersite) then
121 mixer%occsize = mixer%occsize + 2*this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets*this%maxneighbors
122 end if
123 end if
124
125 if (states_are_real(st)) then
126 safe_allocate(mixer%dtmp_occ(1:mixer%occsize, 1))
127 else
128 safe_allocate(mixer%ztmp_occ(1:mixer%occsize, 1))
129 end if
130
131 if (this%level == dft_u_acbn0) then
132 safe_allocate(mixer%tmpU(1:this%norbsets, 1))
133 end if
134
135 pop_sub(lda_u_mixer_init)
136 end subroutine lda_u_mixer_init
137
138 ! ---------------------------------------------------------
139 subroutine lda_u_mixer_clear(mixer, smix)
140 type(lda_u_mixer_t), intent(inout) :: mixer
141 type(mix_t), intent(inout) :: smix
143 if (.not. mixer%apply) return
146 call mixfield_clear(mix_scheme(smix), mixer%mixfield_occ)
147 if (mixer%mixU) call mixfield_clear(mix_scheme(smix), mixer%mixfield_U)
149 pop_sub(lda_u_mixer_clear)
150 end subroutine lda_u_mixer_clear
151
152 ! ---------------------------------------------------------
153 subroutine lda_u_mixer_end(mixer, smix)
154 type(lda_u_mixer_t), intent(inout) :: mixer
155 type(mix_t), intent(inout) :: smix
157 if (.not. mixer%apply) return
158 push_sub(lda_u_mixer_end)
159
160 safe_deallocate_a(mixer%dtmp_occ)
161 safe_deallocate_a(mixer%ztmp_occ)
162 safe_deallocate_a(mixer%tmpU)
163
164 call mixfield_end(smix,mixer%mixfield_occ)
165 if (mixer%mixU) call mixfield_end(smix, mixer%mixfield_U)
166
167 pop_sub(lda_u_mixer_end)
168 end subroutine lda_u_mixer_end
169
170 ! ---------------------------------------------------------
171 subroutine lda_u_mixer_set_vout(this, mixer)
172 type(lda_u_t), intent(in) :: this
173 type(lda_u_mixer_t), intent(inout) :: mixer
174
175 if (this%level == dft_u_none) return
176 push_sub(lda_u_mixer_set_vout)
177
178 if (mixer%realstates) then
179 call dlda_u_get_occupations(this, mixer%dtmp_occ(1:mixer%occsize, 1))
180 call mixfield_set_vout(mixer%mixfield_occ, mixer%dtmp_occ)
181 else
182 call zlda_u_get_occupations(this, mixer%ztmp_occ(1:mixer%occsize, 1))
183 call mixfield_set_vout(mixer%mixfield_occ, mixer%ztmp_occ)
184 end if
185
186 if (this%level == dft_u_acbn0) then
187 call lda_u_get_effectiveu(this, mixer%tmpU(1:this%norbsets, 1))
188 call mixfield_set_vout(mixer%mixfield_U, mixer%tmpU)
189 end if
190
191 pop_sub(lda_u_mixer_set_vout)
192 end subroutine lda_u_mixer_set_vout
193
194 ! ---------------------------------------------------------
195 subroutine lda_u_mixer_set_vin(this, mixer)
196 type(lda_u_t), intent(in) :: this
197 type(lda_u_mixer_t), intent(inout) :: mixer
198
199 if (.not. mixer%apply) return
201
202 if (mixer%realstates) then
203 call dlda_u_get_occupations(this, mixer%dtmp_occ(1:mixer%occsize, 1))
204 call mixfield_set_vin(mixer%mixfield_occ, mixer%dtmp_occ)
205 else
206 call zlda_u_get_occupations(this, mixer%ztmp_occ(1:mixer%occsize, 1))
207 call mixfield_set_vin(mixer%mixfield_occ, mixer%ztmp_occ)
208 end if
209
210 if (this%level == dft_u_acbn0) then
211 call lda_u_get_effectiveu(this, mixer%tmpU(1:this%norbsets, 1))
212 call mixfield_set_vin(mixer%mixfield_U, mixer%tmpU)
213 end if
214
215 pop_sub(lda_u_mixer_set_vin)
216 end subroutine lda_u_mixer_set_vin
217
218 ! ---------------------------------------------------------
219 subroutine lda_u_mixer_get_vnew(this, mixer, st)
220 type(lda_u_t), intent(inout) :: this
221 type(lda_u_mixer_t), intent(inout) :: mixer
222 type(states_elec_t), intent(in) :: st
223
224 if (.not. mixer%apply) return
225 push_sub(lda_u_mixer_get_vnew)
226
227 if (this%level == dft_u_acbn0) then
228 call mixfield_get_vnew(mixer%mixfield_U, mixer%tmpU)
229 call lda_u_set_effectiveu(this, mixer%tmpU(1:this%norbsets, 1))
230 end if
231
233 if (mixer%realstates) then
234 call mixfield_get_vnew(mixer%mixfield_occ, mixer%dtmp_occ)
235 call dlda_u_set_occupations(this, mixer%dtmp_occ(1:mixer%occsize, 1))
236 call dlda_u_update_potential(this, st)
237 else
238 call mixfield_get_vnew(mixer%mixfield_occ, mixer%ztmp_occ)
239 call zlda_u_set_occupations(this, mixer%ztmp_occ(1:mixer%occsize, 1))
240 call zlda_u_update_potential(this, st)
241 end if
242
243 pop_sub(lda_u_mixer_get_vnew)
244 end subroutine lda_u_mixer_get_vnew
245
247end module lda_u_mixer_oct_m
subroutine, public lda_u_mixer_set_vin(this, mixer)
subroutine, public lda_u_mixer_init(this, mixer, st)
subroutine, public lda_u_mixer_clear(mixer, smix)
subroutine, public lda_u_mixer_init_auxmixer(this, namespace, mixer, smix, st)
subroutine, public lda_u_mixer_get_vnew(this, mixer, st)
subroutine, public lda_u_mixer_set_vout(this, mixer)
subroutine, public lda_u_mixer_end(mixer, smix)
subroutine, public lda_u_get_effectiveu(this, Ueff)
Definition: lda_u.F90:916
subroutine, public lda_u_set_effectiveu(this, Ueff)
Definition: lda_u.F90:900
integer, parameter, public dft_u_none
Definition: lda_u.F90:200
subroutine, public dlda_u_get_occupations(this, occ)
Definition: lda_u.F90:3494
subroutine, public zlda_u_update_potential(this, st)
This routine computes the potential that, once multiplied by the projector Pmm' and summed over m and...
Definition: lda_u.F90:4288
subroutine, public dlda_u_set_occupations(this, occ)
Definition: lda_u.F90:3437
integer, parameter, public dft_u_acbn0
Definition: lda_u.F90:200
subroutine, public zlda_u_get_occupations(this, occ)
Definition: lda_u.F90:5532
subroutine, public dlda_u_update_potential(this, st)
This routine computes the potential that, once multiplied by the projector Pmm' and summed over m and...
Definition: lda_u.F90:2301
subroutine, public zlda_u_set_occupations(this, occ)
Definition: lda_u.F90:5475
subroutine, public mixfield_end(smix, mixfield)
Deallocate all arrays of a mixfield instance.
Definition: mix.F90:870
integer pure function, public mix_scheme(this)
Definition: mix.F90:772
subroutine, public mixfield_init(smix, mixfield, d1, d2, d3, func_type)
Initialise all attributes of a mixfield instance.
Definition: mix.F90:826
subroutine, public mix_add_auxmixfield(namespace, smix, mixfield)
Definition: mix.F90:807
subroutine, public mixfield_clear(scheme, mixfield)
Zero all potential and field attributes of a mixfield instance.
Definition: mix.F90:907
integer pure function, public mix_d3(this)
Definition: mix.F90:778
pure logical function, public states_are_real(st)
type(type_t), public type_float
Definition: types.F90:133
type(type_t), public type_cmplx
Definition: types.F90:134
Class to describe DFT+U parameters.
Definition: lda_u.F90:213
God class for mixing.
Definition: mix.F90:213
The states_elec_t class contains all electronic wave functions.
int true(void)