Octopus
lda_u_io.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#include "global.h"
19
20module lda_u_io_oct_m
22 use debug_oct_m
23 use global_oct_m
24 use io_oct_m
26 use ions_oct_m
28 use lda_u_oct_m
29 use mesh_oct_m
32 use mpi_oct_m
36 use parser_oct_m
39 use space_oct_m
43 use unit_oct_m
45
46 implicit none
47
48 private
49
50 public :: &
57 lda_u_load, &
59
60contains
61
63 subroutine lda_u_write_occupation_matrices(dir, this, st, namespace)
64 type(lda_u_t), intent(in) :: this
65 character(len=*), intent(in) :: dir
66 type(states_elec_t), intent(in) :: st
67 type(namespace_t), intent(in) :: namespace
68
69 integer :: iunit, ios, ispin, im, imp, ios2, inn
70 type(orbitalset_t), pointer :: os, os2
71
73
74 if (mpi_grp_is_root(mpi_world)) then ! this the absolute master writes
75 iunit = io_open(trim(dir) // "/occ_matrices", namespace, action='write')
76 write(iunit,'(a)') ' Occupation matrices '
77
78 do ios = 1, this%norbsets
79 os => this%orbsets(ios)
80
81 do ispin = 1,st%d%nspin
82 write(iunit,'(a, i4, a, i4)') ' Orbital set ', ios, ' spin ', ispin
83
84 do im = 1, os%norbs
85 write(iunit,'(1x)',advance='no')
86
87 if (states_are_real(st)) then
88 do imp = 1, os%norbs-1
89 write(iunit,'(f14.8)', advance='no') this%dn(im, imp, ispin, ios)
90 end do
91 write(iunit,'(f14.8)') this%dn(im, os%norbs, ispin, ios)
92 else
93 do imp = 1, os%norbs-1
94 write(iunit,'(f14.8,f14.8)', advance='no') this%zn(im, imp, ispin, ios)
95 end do
96 write(iunit,'(f14.8,f14.8)') this%zn(im, os%norbs, ispin, ios)
97 end if
98 end do
99 end do !ispin
100 end do !iatom
101 call io_close(iunit)
102
103 if (this%level == dft_u_acbn0) then
104 iunit = io_open(trim(dir) // "/renorm_occ_matrices", namespace, action='write')
105 write(iunit,'(a)') ' Renormalized occupation matrices '
106
107 do ios = 1, this%norbsets
108 os => this%orbsets(ios)
109 do ispin = 1, st%d%nspin
110 write(iunit,'(a, i4, a, i4)') ' Orbital set ', ios, ' spin ', ispin
111 do im = 1, os%norbs
112
113 write(iunit,'(1x)',advance='no')
114
115 if (states_are_real(st)) then
116 do imp = 1, os%norbs-1
117 write(iunit,'(f14.8)', advance='no') this%dn_alt(im, imp, ispin, ios)
118 end do
119 write(iunit,'(f14.8)') this%dn_alt(im, os%norbs, ispin, ios)
120 else
121 do imp = 1, os%norbs-1
122 write(iunit,'(f14.8,f14.8)', advance='no') this%zn_alt(im, imp, ispin, ios)
123 end do
124 write(iunit,'(f14.8,f14.8)') this%zn_alt(im, os%norbs, ispin, ios)
125 end if
126 end do
127 end do !ispin
128 end do !iatom
129 call io_close(iunit)
130 end if
131
132 if(this%intersite) then
133 iunit = io_open(trim(dir) // "/intersite_occ_matrices", namespace, action='write')
134 write(iunit,'(a)') ' Intersite occupation matrices '
135 do ios = 1, this%norbsets
136 os => this%orbsets(ios)
137 do ispin = 1, st%d%nspin
138 write(iunit,'(a, i4, a, i4)') ' Orbital set ', ios, ' spin ', ispin
139 do inn = 1, os%nneighbors
140 write(iunit,'(a, i4)') ' Neighbour ', inn
141 ios2 = os%map_os(inn)
142 os2 => this%orbsets(ios2)
143 do im = 1, os%norbs
144
145 write(iunit,'(1x)',advance='no')
146
147 if (states_are_real(st)) then
148 do imp = 1, os2%norbs - 1
149 write(iunit,'(f14.8)', advance='no') this%dn_ij(im, imp, ispin, ios, inn)
150 end do
151 write(iunit,'(f14.8)') this%dn_ij(im, os2%norbs, ispin, ios, inn)
152 else
153 do imp = 1, os2%norbs - 1
154 write(iunit,'(f14.8,f14.8)', advance='no') this%zn_ij(im, imp, ispin, ios, inn)
155 end do
156 write(iunit,'(f14.8,f14.8)') this%zn_ij(im, os2%norbs, ispin, ios, inn)
157 end if
158 end do
159 end do !inn
160 end do !ispin
161 end do ! ios
162 call io_close(iunit)
163 end if
164
165 end if
166
169
170 !---------------------------------------------------------
171 subroutine lda_u_write_effectiveu(dir, this, namespace)
172 type(lda_u_t), intent(in) :: this
173 character(len=*), intent(in) :: dir
174 type(namespace_t), intent(in) :: namespace
175
176 integer :: iunit, ios
177
178 push_sub(lda_u_write_effectiveu)
179
180 if (mpi_grp_is_root(mpi_world)) then ! this the absolute master writes
181 iunit = io_open(trim(dir) // "/effectiveU", namespace, action='write')
182 call lda_u_write_u(this, iunit)
183
184 write(iunit, '(a,a,a)') 'Hubbard U [', trim(units_abbrev(units_out%energy)),']:'
185 write(iunit,'(a,6x,14x,a)') ' Orbital', 'U'
186 do ios = 1, this%norbsets
187 if (.not. this%basisfromstates) then
188 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals) then
189 if (this%orbsets(ios)%nn /= 0) then
190 write(iunit,'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
191 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
192 units_from_atomic(units_out%energy, this%orbsets(ios)%Ubar)
193 else
194 write(iunit,'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
195 l_notation(this%orbsets(ios)%ll), &
196 units_from_atomic(units_out%energy, this%orbsets(ios)%Ubar)
197 end if
198 else
199 if (this%orbsets(ios)%nn /= 0) then
200 write(iunit,'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
201 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
202 int(m_two*(this%orbsets(ios)%jj)), '/2', &
203 units_from_atomic(units_out%energy, this%orbsets(ios)%Ubar)
204 else
205 write(iunit,'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
206 l_notation(this%orbsets(ios)%ll), &
207 int(m_two*(this%orbsets(ios)%jj)), '/2', &
208 units_from_atomic(units_out%energy, this%orbsets(ios)%Ubar)
209 end if
210 end if
211 else
212 write(iunit,'(i4,a10, 3x, f15.6)') ios, 'states', units_from_atomic(units_out%energy, this%orbsets(ios)%Ubar)
213 end if
214 end do
215
216
217 write(iunit, '(a,a,a)') 'Hund J [', trim(units_abbrev(units_out%energy)),']:'
218 write(iunit,'(a,6x,14x,a)') ' Orbital', 'J'
219 do ios = 1, this%norbsets
220 if (.not. this%basisfromstates) then
221 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals) then
222 if (this%orbsets(ios)%nn /= 0) then
223 write(iunit,'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
224 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
225 units_from_atomic(units_out%energy, this%orbsets(ios)%Jbar)
226 else
227 write(iunit,'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
228 l_notation(this%orbsets(ios)%ll), &
229 units_from_atomic(units_out%energy, this%orbsets(ios)%Jbar)
230 end if
231 else
232 if (this%orbsets(ios)%nn /= 0) then
233 write(iunit,'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
234 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
235 int(m_two*(this%orbsets(ios)%jj)), '/2', &
236 units_from_atomic(units_out%energy, this%orbsets(ios)%Jbar)
237 else
238 write(iunit,'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
239 l_notation(this%orbsets(ios)%ll), &
240 int(m_two*(this%orbsets(ios)%jj)), '/2', &
241 units_from_atomic(units_out%energy, this%orbsets(ios)%Jbar)
242 end if
243 end if
244 else
245 write(iunit,'(i4,a10, f15.6)') ios, 'states', units_from_atomic(units_out%energy, this%orbsets(ios)%Jbar)
246 end if
247 end do
248
249 call io_close(iunit)
250 end if
251
253 end subroutine lda_u_write_effectiveu
254
255 !---------------------------------------------------------
256 subroutine lda_u_write_kanamoriu(dir, st, this, namespace)
257 type(lda_u_t), intent(in) :: this
258 type(states_elec_t), intent(in) :: st
259 character(len=*), intent(in) :: dir
260 type(namespace_t), intent(in) :: namespace
261
262 integer :: iunit, ios
263 real(real64), allocatable :: kanamori(:,:)
265 push_sub(lda_u_write_kanamoriu)
266
267 if (mpi_grp_is_root(mpi_world)) then ! this the absolute master writes
268 safe_allocate(kanamori(1:3,1:this%norbsets))
269
270 call compute_acbno_u_kanamori(this, st, kanamori)
271
272 iunit = io_open(trim(dir) // "/kanamoriU", namespace, action='write')
273
274 write(iunit, '(a,a,a)') 'Intraorbital U [', trim(units_abbrev(units_out%energy)),']:'
275 write(iunit,'(a,6x,14x,a)') ' Orbital', 'U'
276 do ios = 1, this%norbsets
277 if (.not. this%basisfromstates) then
278 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals) then
279 if (this%orbsets(ios)%nn /= 0) then
280 write(iunit,'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
281 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
282 units_from_atomic(units_out%energy, kanamori(1,ios))
283 else
284 write(iunit,'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
285 l_notation(this%orbsets(ios)%ll), &
286 units_from_atomic(units_out%energy, kanamori(1,ios))
287 end if
288 else
289 if (this%orbsets(ios)%nn /= 0) then
290 write(iunit,'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
291 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
292 int(m_two*(this%orbsets(ios)%jj)), '/2', &
293 units_from_atomic(units_out%energy, kanamori(1,ios))
294 else
295 write(iunit,'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
296 l_notation(this%orbsets(ios)%ll), &
297 int(m_two*(this%orbsets(ios)%jj)), '/2', &
298 units_from_atomic(units_out%energy, kanamori(1,ios))
299 end if
300 end if
301 else
302 write(iunit,'(i4,a10, 3x, f15.6)') ios, 'states', units_from_atomic(units_out%energy, kanamori(1,ios))
303 end if
304 end do
305
306
307 write(iunit, '(a,a,a)') 'Interorbital Up [', trim(units_abbrev(units_out%energy)),']:'
308 write(iunit,'(a,6x,14x,a)') ' Orbital', 'Up'
309 do ios = 1, this%norbsets
310 if (.not. this%basisfromstates) then
311 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals) then
312 if (this%orbsets(ios)%nn /= 0) then
313 write(iunit,'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
314 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
315 units_from_atomic(units_out%energy, kanamori(2,ios))
316 else
317 write(iunit,'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
318 l_notation(this%orbsets(ios)%ll), &
319 units_from_atomic(units_out%energy, kanamori(2,ios))
320 end if
321 else
322 if (this%orbsets(ios)%nn /= 0) then
323 write(iunit,'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
324 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
325 int(m_two*(this%orbsets(ios)%jj)), '/2', &
326 units_from_atomic(units_out%energy, kanamori(2,ios))
327 else
328 write(iunit,'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
329 l_notation(this%orbsets(ios)%ll), &
330 int(m_two*(this%orbsets(ios)%jj)), '/2', &
331 units_from_atomic(units_out%energy, kanamori(2,ios))
332 end if
333 end if
334 else
335 write(iunit,'(i4,a10, f15.6)') ios, 'states', units_from_atomic(units_out%energy, kanamori(2,ios))
336 end if
337 end do
338
339 write(iunit, '(a,a,a)') 'Hund J [', trim(units_abbrev(units_out%energy)),']:'
340 write(iunit,'(a,6x,14x,a)') ' Orbital', 'J'
341 do ios = 1, this%norbsets
342 if (.not. this%basisfromstates) then
343 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals) then
344 if (this%orbsets(ios)%nn /= 0) then
345 write(iunit,'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
346 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
347 units_from_atomic(units_out%energy, kanamori(3,ios))
348 else
349 write(iunit,'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
350 l_notation(this%orbsets(ios)%ll), &
351 units_from_atomic(units_out%energy, kanamori(3,ios))
352 end if
353 else
354 if (this%orbsets(ios)%nn /= 0) then
355 write(iunit,'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
356 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
357 int(m_two*(this%orbsets(ios)%jj)), '/2', &
358 units_from_atomic(units_out%energy, kanamori(3,ios))
359 else
360 write(iunit,'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
361 l_notation(this%orbsets(ios)%ll), &
362 int(m_two*(this%orbsets(ios)%jj)), '/2', &
363 units_from_atomic(units_out%energy, kanamori(3,ios))
364 end if
365 end if
366 else
367 write(iunit,'(i4,a10, f15.6)') ios, 'states', units_from_atomic(units_out%energy, kanamori(3,ios))
368 end if
369 end do
370
371
372 call io_close(iunit)
373
374 safe_deallocate_a(kanamori)
375 end if
376
377
378 pop_sub(lda_u_write_kanamoriu)
379 end subroutine lda_u_write_kanamoriu
380
381
382
383 !---------------------------------------------------------
384 subroutine lda_u_write_magnetization(dir, this, ions, mesh, st, namespace)
385 type(lda_u_t), intent(in) :: this
386 character(len=*), intent(in) :: dir
387 type(ions_t), intent(in) :: ions
388 class(mesh_t), intent(in) :: mesh
389 type(states_elec_t), intent(in) :: st
390 type(namespace_t), intent(in) :: namespace
391
392 integer :: iunit, ia, ios, im
393 real(real64), allocatable :: mm(:,:)
394
395 if (.not. mpi_grp_is_root(mpi_world)) return
396
398
399 call io_mkdir(dir, namespace)
400 iunit = io_open(trim(dir)//"/magnetization.xsf", namespace, action='write', position='asis')
401
402 if (this%nspins > 1) then
403 safe_allocate(mm(1:3, 1:ions%natoms))
404 mm = m_zero
405 !We compute the magnetization vector for each orbital set
406 do ios = 1, this%norbsets
407 ia = this%orbsets(ios)%iatom
408 do im = 1, this%orbsets(ios)%norbs
409 if (states_are_real(st)) then
410 mm(3, ia) = mm(3, ia) + this%dn(im,im,1,ios) - this%dn(im,im,2,ios)
411 else
412 mm(3, ia) = mm(3, ia) + real(this%zn(im,im,1,ios) - this%zn(im,im,2,ios), real64)
413 !Spinors
414 if (this%nspins /= this%spin_channels) then
415 mm(1, ia) = mm(1, ia) + 2*real(this%zn(im,im,3,ios), real64)
416 mm(2, ia) = mm(2, ia) - 2*aimag(this%zn(im,im,3,ios))
417 end if
418 end if
419 end do !im
420 end do ! ios
421 call write_xsf_geometry(iunit, ions%space, ions%latt, ions%pos, ions%atom, mesh, forces = mm)
422 safe_deallocate_a(mm)
423 end if
424
425 call io_close(iunit)
426
428 end subroutine lda_u_write_magnetization
429
430 !---------------------------------------------------------
431 subroutine lda_u_write_u(this, iunit, namespace)
432 type(lda_u_t), intent(in) :: this
433 integer, optional, intent(in) :: iunit
434 type(namespace_t), optional, intent(in) :: namespace
435
436 integer :: ios
437
438 push_sub(lda_u_write_u)
439
440 write(message(1), '(a,a,a)') 'Effective Hubbard U [', trim(units_abbrev(units_out%energy)),']:'
441 write(message(2),'(a,6x,14x,a)') ' Orbital', 'U'
442 call messages_info(2, iunit=iunit, namespace=namespace)
443
444 do ios = 1, this%norbsets
445 if (.not. this%basisfromstates) then
446 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals) then
447 if (this%orbsets(ios)%nn /= 0) then
448 write(message(1),'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
449 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
450 units_from_atomic(units_out%energy, this%orbsets(ios)%Ueff)
451 else
452 write(message(1),'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
453 l_notation(this%orbsets(ios)%ll), &
454 units_from_atomic(units_out%energy, this%orbsets(ios)%Ueff)
455 end if
456 else
457 if (this%orbsets(ios)%nn /= 0) then
458 write(message(1),'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
459 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
460 int(m_two*(this%orbsets(ios)%jj)), '/2', &
461 units_from_atomic(units_out%energy, this%orbsets(ios)%Ueff)
462 else
463 write(message(1),'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
464 l_notation(this%orbsets(ios)%ll), &
465 int(m_two*(this%orbsets(ios)%jj)), '/2', &
466 units_from_atomic(units_out%energy, this%orbsets(ios)%Ueff)
467 end if
468 end if
469 else
470 write(message(1),'(i4,a10, f15.6)') ios, 'states', units_from_atomic(units_out%energy, this%orbsets(ios)%Ueff)
471 end if
472 call messages_info(1, iunit=iunit, namespace=namespace)
473 end do
474
475 pop_sub(lda_u_write_u)
476 end subroutine lda_u_write_u
478 !---------------------------------------------------------
479 subroutine lda_u_write_v(this, iunit, namespace)
480 type(lda_u_t), intent(in) :: this
481 integer, optional, intent(in) :: iunit
482 type(namespace_t), optional, intent(in) :: namespace
483
484 integer :: ios, icopies, ios2
485
486 if (.not. this%intersite) return
487
488 push_sub(lda_u_write_v)
489
490 write(message(1), '(a,a,a)') 'Effective intersite V [', trim(units_abbrev(units_out%energy)),']:'
491 write(message(2),'(a,14x,a)') ' Orbital', 'V'
492 call messages_info(2, iunit=iunit, namespace=namespace)
493
494 do ios = 1, this%norbsets
495 do icopies = 1, this%orbsets(ios)%nneighbors
496 ios2 = this%orbsets(ios)%map_os(icopies)
497 if (.not.this%basisfromstates) then
498 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals) then
499 if (this%orbsets(ios)%nn /= 0) then
500 write(message(1),'(i4,a10, 2x, i1, a1, i4, 1x, i1, a1, f7.3, f15.6)') ios, &
501 trim(this%orbsets(ios)%spec%get_label()), this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), ios2, &
502 this%orbsets(ios2)%nn, l_notation(this%orbsets(ios2)%ll), &
503 units_from_atomic(units_out%length, this%orbsets(ios)%V_ij(icopies,3+1)), &
504 units_from_atomic(units_out%energy, this%orbsets(ios)%V_ij(icopies,0))
505 call messages_info(1, iunit=iunit, namespace=namespace)
506 else
507 write(message(1),'(i4,a10, 3x, a1, i4, 1x, a1, f7.3, f15.6)') ios, &
508 trim(this%orbsets(ios)%spec%get_label()), l_notation(this%orbsets(ios)%ll), ios2, &
509 l_notation(this%orbsets(ios2)%ll), units_from_atomic(units_out%length, this%orbsets(ios)%V_ij(icopies,3+1)), &
510 units_from_atomic(units_out%energy, this%orbsets(ios)%V_ij(icopies,0))
511 call messages_info(1, iunit=iunit, namespace=namespace)
512 end if
513 else
514 if (this%orbsets(ios)%nn /= 0) then
515 write(message(1),'(i4,a10, 2x, i1, a1, i1, a2, i4, f7.3, f15.6)') ios, &
516 trim(this%orbsets(ios)%spec%get_label()), this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
517 int(m_two*(this%orbsets(ios)%jj)), '/2', ios2, &
518 units_from_atomic(units_out%length, this%orbsets(ios)%V_ij(icopies,3+1)), &
519 units_from_atomic(units_out%energy, this%orbsets(ios)%V_ij(icopies,0))
520 call messages_info(1, iunit=iunit, namespace=namespace)
521 else
522 write(message(1),'(i4,a10, 3x, a1, i1, a2, i4, f7.3, f15.6)') ios, &
523 trim(this%orbsets(ios)%spec%get_label()), l_notation(this%orbsets(ios)%ll), &
524 int(m_two*(this%orbsets(ios)%jj)), '/2', ios2, &
525 units_from_atomic(units_out%length, this%orbsets(ios)%V_ij(icopies,3+1)), &
526 units_from_atomic(units_out%energy, this%orbsets(ios)%V_ij(icopies,0))
527 call messages_info(1, iunit=iunit, namespace=namespace)
528 end if ! this%orbsets(ios)%nn /= 0
529 end if ! this%orbsets(ios)%ndim == 1
530 else
531 write(message(1),'(i4,a10, i4, f7.3, f15.6)') ios, 'states', ios2, &
532 units_from_atomic(units_out%length, this%orbsets(ios)%V_ij(icopies,3+1)), &
533 units_from_atomic(units_out%energy, this%orbsets(ios)%V_ij(icopies,0))
534 call messages_info(1, iunit=iunit, namespace=namespace)
535 end if
536
537 end do ! icopies
538 end do ! ios
539
540 pop_sub(lda_u_write_v)
541 end subroutine lda_u_write_v
542
543
544 ! ---------------------------------------------------------
545 subroutine lda_u_dump(restart, namespace, this, st, mesh, ierr)
546 type(restart_t), intent(in) :: restart
547 type(namespace_t), intent(in) :: namespace
548 type(lda_u_t), intent(in) :: this
549 type(states_elec_t), intent(in) :: st
550 class(mesh_t), intent(in) :: mesh
551 integer, intent(out) :: ierr
552
553 integer :: err, occsize, ios, ncount
554 real(real64), allocatable :: ueff(:), docc(:), veff(:)
555 complex(real64), allocatable :: zocc(:)
556
557 push_sub(lda_u_dump)
558
559 ierr = 0
560
561 if (restart_skip(restart)) then
562 pop_sub(lda_u_dump)
563 return
564 end if
565
566 message(1) = "Debug: Writing DFT+U restart."
567 call messages_info(1, namespace=namespace, debug_only=.true.)
568
569 occsize = this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets
570 if (this%level == dft_u_acbn0) then
571 occsize = occsize*2
572 if (this%intersite) then
573 occsize = occsize + 2*this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets*this%maxneighbors
574 end if
575 end if
576
577
578 if (states_are_real(st)) then
579 safe_allocate(docc(1:occsize))
580 docc = m_zero
581 call dlda_u_get_occupations(this, docc)
582 call drestart_write_binary(restart, "lda_u_occ", occsize, docc, err)
583 if (err /= 0) ierr = ierr + 1
584 safe_deallocate_a(docc)
585 else
586 safe_allocate(zocc(1:occsize))
587 zocc = m_zero
588 call zlda_u_get_occupations(this, zocc)
589 call zrestart_write_binary(restart, "lda_u_occ", occsize, zocc, err)
590 if (err /= 0) ierr = ierr + 1
591 safe_deallocate_a(zocc)
592 end if
593
594
595 if (this%level == dft_u_acbn0) then
596 safe_allocate(ueff(1:this%norbsets))
597 ueff = m_zero
598 call lda_u_get_effectiveu(this, ueff(:))
599 call drestart_write_binary(restart, "lda_u_Ueff", this%norbsets, ueff, err)
600 safe_deallocate_a(ueff)
601 if (err /= 0) ierr = ierr + 1
602
603 if (this%intersite .and. this%maxneighbors > 0) then
604 ncount = 0
605 do ios = 1, this%norbsets
606 ncount = ncount + this%orbsets(ios)%nneighbors
607 end do
608 safe_allocate(veff(1:ncount))
609 call lda_u_get_effectivev(this, veff(:))
610 call drestart_write_binary(restart, "lda_u_Veff", ncount, veff, err)
611 safe_deallocate_a(veff)
612 if (err /= 0) ierr = ierr + 1
613 end if
614 end if
615
616 call lda_u_dump_coulomb_integrals(this, namespace, restart, st, mesh, ierr)
617
618 message(1) = "Debug: Writing DFT+U restart done."
619 call messages_info(1, namespace=namespace, debug_only=.true.)
620
621 pop_sub(lda_u_dump)
622 end subroutine lda_u_dump
623
624
625 ! ---------------------------------------------------------
626 subroutine lda_u_load(restart, this, st, dftu_energy, ierr, occ_only, u_only)
627 type(restart_t), intent(in) :: restart
628 type(lda_u_t), intent(inout) :: this
629 type(states_elec_t), intent(in) :: st
630 real(real64), intent(out) :: dftu_energy
631 integer, intent(out) :: ierr
632 logical, optional, intent(in) :: occ_only
633 logical, optional, intent(in) :: u_only
634
635 integer :: err, occsize, ncount, ios
636 real(real64), allocatable :: ueff(:), docc(:), veff(:)
637 complex(real64), allocatable :: zocc(:)
639 push_sub(lda_u_load)
640
641 ierr = 0
642
643 if (restart_skip(restart)) then
644 pop_sub(lda_u_load)
645 return
646 end if
647
648 message(1) = "Debug: Reading DFT+U restart."
649 call messages_info(1, debug_only=.true.)
650
651 !We have to read the effective U first, as we call lda_u_uptade_potential latter
652 if (this%level == dft_u_acbn0 .and. .not. optional_default(occ_only, .false.)) then
653 safe_allocate(ueff(1:this%norbsets))
654 call drestart_read_binary(restart, "lda_u_Ueff", this%norbsets, ueff, err)
655 if (err /= 0) then
656 ierr = ierr + 1
657 ueff = m_zero
658 end if
659 call lda_u_set_effectiveu(this, ueff)
660 safe_deallocate_a(ueff)
661
662 if (this%intersite .and. this%maxneighbors > 0) then
663 ncount = 0
664 do ios = 1, this%norbsets
665 ncount = ncount + this%orbsets(ios)%nneighbors
666 end do
667 safe_allocate(veff(1:ncount))
668 call drestart_read_binary(restart, "lda_u_Veff", ncount, veff, err)
669 if (err /= 0) then
670 ierr = ierr + 1
671 veff = m_zero
672 end if
673 call lda_u_set_effectivev(this, veff)
674 safe_deallocate_a(veff)
675 end if
676 end if
677
678
679 if (.not. optional_default(u_only, .false.)) then
680 occsize = this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets
681 if (this%level == dft_u_acbn0) then
682 occsize = occsize*2
683 if (this%intersite) then
684 occsize = occsize + 2*this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets*this%maxneighbors
685 end if
686 end if
687
688
689 if (states_are_real(st)) then
690 safe_allocate(docc(1:occsize))
691 call drestart_read_binary(restart, "lda_u_occ", occsize, docc, err)
692 if (err /= 0) then
693 ierr = ierr + 1
694 docc = m_zero
695 end if
696 call dlda_u_set_occupations(this, docc)
697 safe_deallocate_a(docc)
698 else
699 safe_allocate(zocc(1:occsize))
700 call zrestart_read_binary(restart, "lda_u_occ", occsize, zocc, err)
701 if (err /= 0) then
702 ierr = ierr + 1
703 zocc = m_z0
704 end if
705 call zlda_u_set_occupations(this, zocc)
706 safe_deallocate_a(zocc)
707 end if
708 end if
709
710 if (states_are_real(st)) then
711 call dcompute_dftu_energy(this, dftu_energy, st)
712 call dlda_u_update_potential(this, st)
713 else
714 call zcompute_dftu_energy(this, dftu_energy, st)
715 call zlda_u_update_potential(this, st)
716 end if
717
718 message(1) = "Debug: Reading DFT+U restart done."
719 call messages_info(1, debug_only=.true.)
720
721 pop_sub(lda_u_load)
722 end subroutine lda_u_load
723
724 ! ---------------------------------------------------------
725 subroutine lda_u_dump_coulomb_integrals(this, namespace, restart, st, mesh, ierr)
726 type(lda_u_t), intent(in) :: this
727 type(namespace_t), intent(in) :: namespace
728 type(restart_t), intent(in) :: restart
729 type(states_elec_t), intent(in) :: st
730 class(mesh_t), intent(in) :: mesh
731 integer, intent(out) :: ierr
732
733 integer :: coulomb_int_file, idim1, idim2, err
734 integer :: ios, im, imp, impp, imppp
735 character(len=256) :: lines(3)
736 logical :: complex_coulomb_integrals
737
739
740 ierr = 0
741
742 if (restart_skip(restart) .or. this%level == dft_u_empirical) then
744 return
745 end if
746
747 message(1) = "Debug: Writing Coulomb integrals restart."
748 call messages_info(1, debug_only=.true.)
749
750 complex_coulomb_integrals = .false.
751 do ios = 1, this%norbsets
752 if (this%orbsets(ios)%ndim > 1) complex_coulomb_integrals = .true.
753 end do
754
755 coulomb_int_file = restart_open(restart, 'coulomb_integrals')
756 ! sanity checks. Example file 'coulomb_int_file':
757 ! norb= 2
758 ! dim= 1
759 ! checksum= xxxxxxxxxxx
760 write(lines(1), '(a20,i21)') "norb=", this%norbsets
761 write(lines(2), '(a20,i21)') "dim=", st%d%dim
762 write(lines(3), '(a20,i21)') "checksum=", mesh%idx%checksum
763 call restart_write(restart, coulomb_int_file, lines, 3, err)
764 if (err /= 0) ierr = ierr - 2
765
766 ios_cycle: do ios = 1, this%norbsets
767 do im = 1, this%orbsets(ios)%norbs
768 do imp = 1, this%orbsets(ios)%norbs
769 do impp = 1, this%orbsets(ios)%norbs
770 do imppp = 1, this%orbsets(ios)%norbs
771 if(.not. complex_coulomb_integrals) then
772 write(lines(1), '(i4,i4,i4,i4,i4,e20.12)') ios, im, imp, impp, imppp, this%coulomb(im, imp, impp, imppp, ios)
773 else
774 do idim1 = 1, st%d%dim
775 do idim2 = 1, st%d%dim
776 write(lines(1), '(i4,i4,i4,i4,i4,i4,i4,2e20.12)') ios, im, imp, impp, imppp, idim1, idim2, &
777 real(this%zcoulomb(im, imp, impp, imppp, idim1, idim2, ios), real64) , &
778 aimag(this%zcoulomb(im, imp, impp, imppp, idim1, idim2, ios))
779 end do
780 end do
781 end if
782 call restart_write(restart, coulomb_int_file, lines, 1, err)
783 if (err /=0) then
784 ierr = ierr - 2**2
785 exit ios_cycle
786 end if
787 end do
788 end do
789 end do
790 end do
791 end do ios_cycle
792 call restart_close(restart, coulomb_int_file)
793
794 message(1) = "Debug: Writing Coulomb integrals restart done."
795 call messages_info(1, debug_only=.true.)
796
798 end subroutine lda_u_dump_coulomb_integrals
799
800end module lda_u_io_oct_m
character(len=1), dimension(0:3), parameter, public l_notation
real(real64), parameter, public m_two
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
complex(real64), parameter, public m_z0
Definition: global.F90:198
subroutine, public write_xsf_geometry(iunit, space, latt, pos, atoms, mesh, forces, index)
for format specification see: http:
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:311
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
subroutine, public lda_u_write_occupation_matrices(dir, this, st, namespace)
Prints the occupation matrices at the end of the scf calculation.
Definition: lda_u_io.F90:157
subroutine, public lda_u_write_kanamoriu(dir, st, this, namespace)
Definition: lda_u_io.F90:350
subroutine, public lda_u_dump(restart, namespace, this, st, mesh, ierr)
Definition: lda_u_io.F90:639
subroutine, public lda_u_write_u(this, iunit, namespace)
Definition: lda_u_io.F90:525
subroutine, public lda_u_load(restart, this, st, dftu_energy, ierr, occ_only, u_only)
Definition: lda_u_io.F90:720
subroutine lda_u_dump_coulomb_integrals(this, namespace, restart, st, mesh, ierr)
Definition: lda_u_io.F90:819
subroutine, public lda_u_write_magnetization(dir, this, ions, mesh, st, namespace)
Definition: lda_u_io.F90:478
subroutine, public lda_u_write_effectiveu(dir, this, namespace)
Definition: lda_u_io.F90:265
subroutine, public lda_u_write_v(this, iunit, namespace)
Definition: lda_u_io.F90:573
subroutine, public lda_u_get_effectiveu(this, Ueff)
Definition: lda_u.F90:917
subroutine, public lda_u_set_effectiveu(this, Ueff)
Definition: lda_u.F90:901
subroutine, public dlda_u_get_occupations(this, occ)
Definition: lda_u.F90:3470
subroutine, public zcompute_dftu_energy(this, energy, st)
This routine computes the value of the double counting term in the DFT+U energy.
Definition: lda_u.F90:4176
integer, parameter, public dft_u_empirical
Definition: lda_u.F90:201
subroutine, public compute_acbno_u_kanamori(this, st, kanamori)
Definition: lda_u.F90:864
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:4250
subroutine, public dlda_u_set_occupations(this, occ)
Definition: lda_u.F90:3413
subroutine, public lda_u_get_effectivev(this, Veff)
Definition: lda_u.F90:951
subroutine, public lda_u_set_effectivev(this, Veff)
Definition: lda_u.F90:933
integer, parameter, public dft_u_acbn0
Definition: lda_u.F90:201
subroutine, public zlda_u_get_occupations(this, occ)
Definition: lda_u.F90:5492
subroutine, public dcompute_dftu_energy(this, energy, st)
This routine computes the value of the double counting term in the DFT+U energy.
Definition: lda_u.F90:2205
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:2279
subroutine, public zlda_u_set_occupations(this, occ)
Definition: lda_u.F90:5435
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
subroutine, public restart_close(restart, iunit)
Close a file previously opened with restart_open.
Definition: restart.F90:952
subroutine, public restart_write(restart, iunit, lines, nlines, ierr)
Definition: restart.F90:908
logical pure function, public restart_skip(restart)
Returns true if the restart information should neither be read nor written. This might happen because...
Definition: restart.F90:969
integer function, public restart_open(restart, filename, status, position, silent)
Open file "filename" found inside the current restart directory. Depending on the type of restart,...
Definition: restart.F90:861
pure logical function, public states_are_real(st)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
Class to describe DFT+U parameters.
Definition: lda_u.F90:214
Describes mesh distribution to nodes.
Definition: mesh.F90:186
The states_elec_t class contains all electronic wave functions.
int true(void)