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,f7.3,a)') 'Hubbard U [', &
185 trim(units_abbrev(units_out%energy)),']:'
186 write(iunit,'(a,6x,14x,a)') ' Orbital', 'U'
187 do ios = 1, this%norbsets
188 if (.not. this%basisfromstates) then
189 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals) then
190 if (this%orbsets(ios)%nn /= 0) then
191 write(iunit,'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
192 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
193 units_from_atomic(units_out%energy, this%orbsets(ios)%Ubar)
194 else
195 write(iunit,'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
196 l_notation(this%orbsets(ios)%ll), &
197 units_from_atomic(units_out%energy, this%orbsets(ios)%Ubar)
198 end if
199 else
200 if (this%orbsets(ios)%nn /= 0) then
201 write(iunit,'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
202 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
203 int(m_two*(this%orbsets(ios)%jj)), '/2', &
204 units_from_atomic(units_out%energy, this%orbsets(ios)%Ubar)
205 else
206 write(iunit,'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
207 l_notation(this%orbsets(ios)%ll), &
208 int(m_two*(this%orbsets(ios)%jj)), '/2', &
209 units_from_atomic(units_out%energy, this%orbsets(ios)%Ubar)
210 end if
211 end if
212 else
213 write(iunit,'(i4,a10, 3x, f15.6)') ios, 'states', units_from_atomic(units_out%energy, this%orbsets(ios)%Ubar)
214 end if
215 end do
216
217
218 write(iunit, '(a,a,a,f7.3,a)') 'Hund J [', &
219 trim(units_abbrev(units_out%energy)),']:'
220 write(iunit,'(a,6x,14x,a)') ' Orbital', 'J'
221 do ios = 1, this%norbsets
222 if (.not. this%basisfromstates) then
223 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals) then
224 if (this%orbsets(ios)%nn /= 0) then
225 write(iunit,'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
226 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
227 units_from_atomic(units_out%energy, this%orbsets(ios)%Jbar)
228 else
229 write(iunit,'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
230 l_notation(this%orbsets(ios)%ll), &
231 units_from_atomic(units_out%energy, this%orbsets(ios)%Jbar)
232 end if
233 else
234 if (this%orbsets(ios)%nn /= 0) then
235 write(iunit,'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
236 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
237 int(m_two*(this%orbsets(ios)%jj)), '/2', &
238 units_from_atomic(units_out%energy, this%orbsets(ios)%Jbar)
239 else
240 write(iunit,'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
241 l_notation(this%orbsets(ios)%ll), &
242 int(m_two*(this%orbsets(ios)%jj)), '/2', &
243 units_from_atomic(units_out%energy, this%orbsets(ios)%Jbar)
244 end if
245 end if
246 else
247 write(iunit,'(i4,a10, f15.6)') ios, 'states', units_from_atomic(units_out%energy, this%orbsets(ios)%Jbar)
248 end if
249 end do
250
251 call io_close(iunit)
252 end if
253
255 end subroutine lda_u_write_effectiveu
256
257 !---------------------------------------------------------
258 subroutine lda_u_write_kanamoriu(dir, st, this, namespace)
259 type(lda_u_t), intent(in) :: this
260 type(states_elec_t), intent(in) :: st
261 character(len=*), intent(in) :: dir
262 type(namespace_t), intent(in) :: namespace
263
264 integer :: iunit, ios
265 real(real64), allocatable :: kanamori(:,:)
266
267 push_sub(lda_u_write_kanamoriu)
268
269 if (mpi_grp_is_root(mpi_world)) then ! this the absolute master writes
270 safe_allocate(kanamori(1:3,1:this%norbsets))
271
272 call compute_acbno_u_kanamori(this, st, kanamori)
273
274 iunit = io_open(trim(dir) // "/kanamoriU", namespace, action='write')
275
276 write(iunit, '(a,a,a,f7.3,a)') 'Intraorbital U [', &
277 trim(units_abbrev(units_out%energy)),']:'
278 write(iunit,'(a,6x,14x,a)') ' Orbital', 'U'
279 do ios = 1, this%norbsets
280 if (.not. this%basisfromstates) then
281 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals) then
282 if (this%orbsets(ios)%nn /= 0) then
283 write(iunit,'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
284 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
285 units_from_atomic(units_out%energy, kanamori(1,ios))
286 else
287 write(iunit,'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
288 l_notation(this%orbsets(ios)%ll), &
289 units_from_atomic(units_out%energy, kanamori(1,ios))
290 end if
291 else
292 if (this%orbsets(ios)%nn /= 0) then
293 write(iunit,'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
294 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
295 int(m_two*(this%orbsets(ios)%jj)), '/2', &
296 units_from_atomic(units_out%energy, kanamori(1,ios))
297 else
298 write(iunit,'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
299 l_notation(this%orbsets(ios)%ll), &
300 int(m_two*(this%orbsets(ios)%jj)), '/2', &
301 units_from_atomic(units_out%energy, kanamori(1,ios))
302 end if
303 end if
304 else
305 write(iunit,'(i4,a10, 3x, f15.6)') ios, 'states', units_from_atomic(units_out%energy, kanamori(1,ios))
306 end if
307 end do
308
309
310 write(iunit, '(a,a,a,f7.3,a)') 'Interorbital Up [', &
311 trim(units_abbrev(units_out%energy)),']:'
312 write(iunit,'(a,6x,14x,a)') ' Orbital', 'Up'
313 do ios = 1, this%norbsets
314 if (.not. this%basisfromstates) then
315 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals) then
316 if (this%orbsets(ios)%nn /= 0) then
317 write(iunit,'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
318 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
319 units_from_atomic(units_out%energy, kanamori(2,ios))
320 else
321 write(iunit,'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
322 l_notation(this%orbsets(ios)%ll), &
323 units_from_atomic(units_out%energy, kanamori(2,ios))
324 end if
325 else
326 if (this%orbsets(ios)%nn /= 0) then
327 write(iunit,'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
328 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
329 int(m_two*(this%orbsets(ios)%jj)), '/2', &
330 units_from_atomic(units_out%energy, kanamori(2,ios))
331 else
332 write(iunit,'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
333 l_notation(this%orbsets(ios)%ll), &
334 int(m_two*(this%orbsets(ios)%jj)), '/2', &
335 units_from_atomic(units_out%energy, kanamori(2,ios))
336 end if
337 end if
338 else
339 write(iunit,'(i4,a10, f15.6)') ios, 'states', units_from_atomic(units_out%energy, kanamori(2,ios))
340 end if
341 end do
342
343 write(iunit, '(a,a,a,f7.3,a)') 'Hund J [', &
344 trim(units_abbrev(units_out%energy)),']:'
345 write(iunit,'(a,6x,14x,a)') ' Orbital', 'J'
346 do ios = 1, this%norbsets
347 if (.not. this%basisfromstates) then
348 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals) then
349 if (this%orbsets(ios)%nn /= 0) then
350 write(iunit,'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
351 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
352 units_from_atomic(units_out%energy, kanamori(3,ios))
353 else
354 write(iunit,'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
355 l_notation(this%orbsets(ios)%ll), &
356 units_from_atomic(units_out%energy, kanamori(3,ios))
357 end if
358 else
359 if (this%orbsets(ios)%nn /= 0) then
360 write(iunit,'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
361 this%orbsets(ios)%nn, 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 else
365 write(iunit,'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
366 l_notation(this%orbsets(ios)%ll), &
367 int(m_two*(this%orbsets(ios)%jj)), '/2', &
368 units_from_atomic(units_out%energy, kanamori(3,ios))
369 end if
370 end if
371 else
372 write(iunit,'(i4,a10, f15.6)') ios, 'states', units_from_atomic(units_out%energy, kanamori(3,ios))
373 end if
374 end do
375
376
377 call io_close(iunit)
378
379 safe_deallocate_a(kanamori)
380 end if
381
382
383 pop_sub(lda_u_write_kanamoriu)
384 end subroutine lda_u_write_kanamoriu
385
386
387
388 !---------------------------------------------------------
389 subroutine lda_u_write_magnetization(dir, this, ions, mesh, st, namespace)
390 type(lda_u_t), intent(in) :: this
391 character(len=*), intent(in) :: dir
392 type(ions_t), intent(in) :: ions
393 class(mesh_t), intent(in) :: mesh
394 type(states_elec_t), intent(in) :: st
395 type(namespace_t), intent(in) :: namespace
396
397 integer :: iunit, ia, ios, im
398 real(real64), allocatable :: mm(:,:)
399
400 if (.not. mpi_grp_is_root(mpi_world)) return
401
403
404 call io_mkdir(dir, namespace)
405 iunit = io_open(trim(dir)//"/magnetization.xsf", namespace, action='write', position='asis')
406
407 if (this%nspins > 1) then
408 safe_allocate(mm(1:mesh%box%dim, 1:ions%natoms))
409 mm = m_zero
410 !We compute the magnetization vector for each orbital set
411 do ios = 1, this%norbsets
412 ia = this%orbsets(ios)%iatom
413 do im = 1, this%orbsets(ios)%norbs
414 if (states_are_real(st)) then
415 mm(3, ia) = mm(3, ia) + this%dn(im,im,1,ios) - this%dn(im,im,2,ios)
416 else
417 mm(3, ia) = mm(3, ia) + real(this%zn(im,im,1,ios) - this%zn(im,im,2,ios), real64)
418 !Spinors
419 if (this%nspins /= this%spin_channels) then
420 mm(1, ia) = mm(1, ia) + 2*real(this%zn(im,im,3,ios), real64)
421 mm(2, ia) = mm(2, ia) - 2*aimag(this%zn(im,im,3,ios))
422 end if
423 end if
424 end do !im
425 end do ! ios
426 call write_xsf_geometry(iunit, ions%space, ions%latt, ions%pos, ions%atom, mesh, forces = mm)
427 safe_deallocate_a(mm)
428 end if
429
430 call io_close(iunit)
431
433 end subroutine lda_u_write_magnetization
434
435 !---------------------------------------------------------
436 subroutine lda_u_write_u(this, iunit, namespace)
437 type(lda_u_t), intent(in) :: this
438 integer, optional, intent(in) :: iunit
439 type(namespace_t), optional, intent(in) :: namespace
440
441 integer :: ios
442
443 push_sub(lda_u_write_u)
444
445 write(message(1), '(a,a,a,f7.3,a)') 'Effective Hubbard U [', &
446 trim(units_abbrev(units_out%energy)),']:'
447 write(message(2),'(a,6x,14x,a)') ' Orbital', 'U'
448 call messages_info(2, iunit=iunit, namespace=namespace)
449
450 do ios = 1, this%norbsets
451 if (.not. this%basisfromstates) then
452 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals) then
453 if (this%orbsets(ios)%nn /= 0) then
454 write(message(1),'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
455 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
456 units_from_atomic(units_out%energy, this%orbsets(ios)%Ueff)
457 else
458 write(message(1),'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
459 l_notation(this%orbsets(ios)%ll), &
460 units_from_atomic(units_out%energy, this%orbsets(ios)%Ueff)
461 end if
462 else
463 if (this%orbsets(ios)%nn /= 0) then
464 write(message(1),'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
465 this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), &
466 int(m_two*(this%orbsets(ios)%jj)), '/2', &
467 units_from_atomic(units_out%energy, this%orbsets(ios)%Ueff)
468 else
469 write(message(1),'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
470 l_notation(this%orbsets(ios)%ll), &
471 int(m_two*(this%orbsets(ios)%jj)), '/2', &
472 units_from_atomic(units_out%energy, this%orbsets(ios)%Ueff)
473 end if
474 end if
475 else
476 write(message(1),'(i4,a10, f15.6)') ios, 'states', units_from_atomic(units_out%energy, this%orbsets(ios)%Ueff)
477 end if
478 call messages_info(1, iunit=iunit, namespace=namespace)
479 end do
480
481 pop_sub(lda_u_write_u)
482 end subroutine lda_u_write_u
483
484 !---------------------------------------------------------
485 subroutine lda_u_write_v(this, iunit, namespace)
486 type(lda_u_t), intent(in) :: this
487 integer, optional, intent(in) :: iunit
488 type(namespace_t), optional, intent(in) :: namespace
489
490 integer :: ios, icopies, ios2
491
492 if (.not. this%intersite) return
493
494 push_sub(lda_u_write_v)
495
496 write(message(1), '(a,a,a,f7.3,a)') 'Effective intersite V [', &
497 trim(units_abbrev(units_out%energy)),']:'
498 write(message(2),'(a,14x,a)') ' Orbital', 'V'
499 call messages_info(2, iunit=iunit, namespace=namespace)
500
501 do ios = 1, this%norbsets
502 do icopies = 1, this%orbsets(ios)%nneighbors
503 ios2 = this%orbsets(ios)%map_os(icopies)
504 if (.not.this%basisfromstates) then
505 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals) then
506 if (this%orbsets(ios)%nn /= 0) then
507 write(message(1),'(i4,a10, 2x, i1, a1, i4, 1x, i1, a1, f7.3, f15.6)') ios, &
508 trim(this%orbsets(ios)%spec%get_label()), this%orbsets(ios)%nn, l_notation(this%orbsets(ios)%ll), ios2, &
509 this%orbsets(ios2)%nn, l_notation(this%orbsets(ios2)%ll), &
510 units_from_atomic(units_out%length, this%orbsets(ios)%V_ij(icopies,3+1)), &
511 units_from_atomic(units_out%energy, this%orbsets(ios)%V_ij(icopies,0))
512 call messages_info(1, iunit=iunit, namespace=namespace)
513 else
514 write(message(1),'(i4,a10, 3x, a1, i4, 1x, a1, f7.3, f15.6)') ios, &
515 trim(this%orbsets(ios)%spec%get_label()), l_notation(this%orbsets(ios)%ll), ios2, &
516 l_notation(this%orbsets(ios2)%ll), units_from_atomic(units_out%length, this%orbsets(ios)%V_ij(icopies,3+1)), &
517 units_from_atomic(units_out%energy, this%orbsets(ios)%V_ij(icopies,0))
518 call messages_info(1, iunit=iunit, namespace=namespace)
519 end if
520 else
521 if (this%orbsets(ios)%nn /= 0) then
522 write(message(1),'(i4,a10, 2x, i1, a1, i1, a2, i4, f7.3, f15.6)') ios, &
523 trim(this%orbsets(ios)%spec%get_label()), this%orbsets(ios)%nn, 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 else
529 write(message(1),'(i4,a10, 3x, a1, i1, a2, i4, f7.3, f15.6)') ios, &
530 trim(this%orbsets(ios)%spec%get_label()), l_notation(this%orbsets(ios)%ll), &
531 int(m_two*(this%orbsets(ios)%jj)), '/2', 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 ! this%orbsets(ios)%nn /= 0
536 end if ! this%orbsets(ios)%ndim == 1
537 else
538 write(message(1),'(i4,a10, i4, f7.3, f15.6)') ios, 'states', ios2, &
539 units_from_atomic(units_out%length, this%orbsets(ios)%V_ij(icopies,3+1)), &
540 units_from_atomic(units_out%energy, this%orbsets(ios)%V_ij(icopies,0))
541 call messages_info(1, iunit=iunit, namespace=namespace)
542 end if
543
544 end do ! icopies
545 end do ! ios
546
547 pop_sub(lda_u_write_v)
548 end subroutine lda_u_write_v
549
550
551 ! ---------------------------------------------------------
552 subroutine lda_u_dump(restart, namespace, this, st, mesh, ierr)
553 type(restart_t), intent(in) :: restart
554 type(namespace_t), intent(in) :: namespace
555 type(lda_u_t), intent(in) :: this
556 type(states_elec_t), intent(in) :: st
557 class(mesh_t), intent(in) :: mesh
558 integer, intent(out) :: ierr
559
560 integer :: err, occsize, ios, ncount
561 real(real64), allocatable :: ueff(:), docc(:), veff(:)
562 complex(real64), allocatable :: zocc(:)
563
564 push_sub(lda_u_dump)
565
566 ierr = 0
567
568 if (restart_skip(restart)) then
569 pop_sub(lda_u_dump)
570 return
571 end if
572
573 if (debug%info) then
574 message(1) = "Debug: Writing DFT+U restart."
575 call messages_info(1, namespace=namespace)
576 end if
577
578 occsize = this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets
579 if (this%level == dft_u_acbn0) then
580 occsize = occsize*2
581 if (this%intersite) then
582 occsize = occsize + 2*this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets*this%maxneighbors
583 end if
584 end if
585
586
587 if (states_are_real(st)) then
588 safe_allocate(docc(1:occsize))
589 docc = m_zero
590 call dlda_u_get_occupations(this, docc)
591 call drestart_write_binary(restart, "lda_u_occ", occsize, docc, err)
592 if (err /= 0) ierr = ierr + 1
593 safe_deallocate_a(docc)
594 else
595 safe_allocate(zocc(1:occsize))
596 zocc = m_zero
597 call zlda_u_get_occupations(this, zocc)
598 call zrestart_write_binary(restart, "lda_u_occ", occsize, zocc, err)
599 if (err /= 0) ierr = ierr + 1
600 safe_deallocate_a(zocc)
601 end if
602
603
604 if (this%level == dft_u_acbn0) then
605 safe_allocate(ueff(1:this%norbsets))
606 ueff = m_zero
607 call lda_u_get_effectiveu(this, ueff(:))
608 call drestart_write_binary(restart, "lda_u_Ueff", this%norbsets, ueff, err)
609 safe_deallocate_a(ueff)
610 if (err /= 0) ierr = ierr + 1
611
612 if (this%intersite .and. this%maxneighbors > 0) then
613 ncount = 0
614 do ios = 1, this%norbsets
615 ncount = ncount + this%orbsets(ios)%nneighbors
616 end do
617 safe_allocate(veff(1:ncount))
618 call lda_u_get_effectivev(this, veff(:))
619 call drestart_write_binary(restart, "lda_u_Veff", ncount, veff, err)
620 safe_deallocate_a(veff)
621 if (err /= 0) ierr = ierr + 1
622 end if
623 end if
624
625 call lda_u_dump_coulomb_integrals(this, namespace, restart, st, mesh, ierr)
626
627 if (debug%info) then
628 message(1) = "Debug: Writing DFT+U restart done."
629 call messages_info(1, namespace=namespace)
630 end if
631
632 pop_sub(lda_u_dump)
633 end subroutine lda_u_dump
634
635
636 ! ---------------------------------------------------------
637 subroutine lda_u_load(restart, this, st, dftu_energy, ierr, occ_only, u_only)
638 type(restart_t), intent(in) :: restart
639 type(lda_u_t), intent(inout) :: this
640 type(states_elec_t), intent(in) :: st
641 real(real64), intent(out) :: dftu_energy
642 integer, intent(out) :: ierr
643 logical, optional, intent(in) :: occ_only
644 logical, optional, intent(in) :: u_only
646 integer :: err, occsize, ncount, ios
647 real(real64), allocatable :: ueff(:), docc(:), veff(:)
648 complex(real64), allocatable :: zocc(:)
649
650 push_sub(lda_u_load)
651
652 ierr = 0
653
654 if (restart_skip(restart)) then
655 ierr = -1
656 pop_sub(lda_u_load)
657 return
658 end if
659
660 if (debug%info) then
661 message(1) = "Debug: Reading DFT+U restart."
662 call messages_info(1)
663 end if
664
665 !We have to read the effective U first, as we call lda_u_uptade_potential latter
666 if (this%level == dft_u_acbn0 .and. .not. optional_default(occ_only, .false.)) then
667 safe_allocate(ueff(1:this%norbsets))
668 call drestart_read_binary(restart, "lda_u_Ueff", this%norbsets, ueff, err)
669 if (err /= 0) then
670 ierr = ierr + 1
671 ueff = m_zero
672 end if
673 call lda_u_set_effectiveu(this, ueff)
674 safe_deallocate_a(ueff)
675
676 if (this%intersite .and. this%maxneighbors > 0) then
677 ncount = 0
678 do ios = 1, this%norbsets
679 ncount = ncount + this%orbsets(ios)%nneighbors
680 end do
681 safe_allocate(veff(1:ncount))
682 call drestart_read_binary(restart, "lda_u_Veff", ncount, veff, err)
683 if (err /= 0) then
684 ierr = ierr + 1
685 veff = m_zero
686 end if
687 call lda_u_set_effectivev(this, veff)
688 safe_deallocate_a(veff)
689 end if
690 end if
691
692
693 if (.not. optional_default(u_only, .false.)) then
694 occsize = this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets
695 if (this%level == dft_u_acbn0) then
696 occsize = occsize*2
697 if (this%intersite) then
698 occsize = occsize + 2*this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets*this%maxneighbors
699 end if
700 end if
701
702
703 if (states_are_real(st)) then
704 safe_allocate(docc(1:occsize))
705 call drestart_read_binary(restart, "lda_u_occ", occsize, docc, err)
706 if (err /= 0) then
707 ierr = ierr + 1
708 docc = m_zero
709 end if
710 call dlda_u_set_occupations(this, docc)
711 safe_deallocate_a(docc)
712 else
713 safe_allocate(zocc(1:occsize))
714 call zrestart_read_binary(restart, "lda_u_occ", occsize, zocc, err)
715 if (err /= 0) then
716 ierr = ierr + 1
717 zocc = m_z0
718 end if
719 call zlda_u_set_occupations(this, zocc)
720 safe_deallocate_a(zocc)
721 end if
722 end if
723
724 if (states_are_real(st)) then
725 call dcompute_dftu_energy(this, dftu_energy, st)
726 call dlda_u_update_potential(this, st)
727 else
728 call zcompute_dftu_energy(this, dftu_energy, st)
729 call zlda_u_update_potential(this, st)
730 end if
731
732 if (debug%info) then
733 message(1) = "Debug: Reading DFT+U restart done."
734 call messages_info(1)
735 end if
736
737 pop_sub(lda_u_load)
738 end subroutine lda_u_load
739
740 ! ---------------------------------------------------------
741 subroutine lda_u_dump_coulomb_integrals(this, namespace, restart, st, mesh, ierr)
742 type(lda_u_t), intent(in) :: this
743 type(namespace_t), intent(in) :: namespace
744 type(restart_t), intent(in) :: restart
745 type(states_elec_t), intent(in) :: st
746 class(mesh_t), intent(in) :: mesh
747 integer, intent(out) :: ierr
748
749 integer :: coulomb_int_file, idim1, idim2, err
750 integer :: ios, im, imp, impp, imppp
751 character(len=256) :: lines(3)
752 logical :: complex_coulomb_integrals
753
755
756 ierr = 0
757
758 if (restart_skip(restart) .or. this%level == dft_u_empirical) then
759 ierr = -1
761 return
762 end if
763
764 if (debug%info) then
765 message(1) = "Debug: Writing Coulomb integrals restart."
766 call messages_info(1)
767 end if
768
769 complex_coulomb_integrals = .false.
770 do ios = 1, this%norbsets
771 if (this%orbsets(ios)%ndim > 1) complex_coulomb_integrals = .true.
772 end do
773
774 coulomb_int_file = restart_open(restart, 'coulomb_integrals')
775 ! sanity checks. Example file 'coulomb_int_file':
776 ! norb= 2
777 ! dim= 1
778 ! checksum= xxxxxxxxxxx
779 write(lines(1), '(a20,i21)') "norb=", this%norbsets
780 write(lines(2), '(a20,i21)') "dim=", st%d%dim
781 write(lines(3), '(a20,i21)') "checksum=", mesh%idx%checksum
782 call restart_write(restart, coulomb_int_file, lines, 3, err)
783 if (err /= 0) ierr = ierr - 2
784
785 do ios = 1, this%norbsets
786 do im = 1, this%orbsets(ios)%norbs
787 do imp = 1, this%orbsets(ios)%norbs
788 do impp = 1, this%orbsets(ios)%norbs
789 do imppp = 1, this%orbsets(ios)%norbs
790 if(.not. complex_coulomb_integrals) then
791 write(lines(1), '(i4,i4,i4,i4,i4,e20.12)') ios, im, imp, impp, imppp, this%coulomb(im, imp, impp, imppp, ios)
792 else
793 do idim1 = 1, st%d%dim
794 do idim2 = 1, st%d%dim
795 write(lines(1), '(i4,i4,i4,i4,i4,i4,i4,2e20.12)') ios, im, imp, impp, imppp, idim1, idim2, &
796 real(this%zcoulomb(im, imp, impp, imppp, idim1, idim2, ios), real64) , &
797 aimag(this%zcoulomb(im, imp, impp, imppp, idim1, idim2, ios))
798 end do
799 end do
800 end if
801 call restart_write(restart, coulomb_int_file, lines, 1, err)
802 if (err /=0) then
803 ierr = ierr - 2**2
804 exit
805 end if
806 end do
807 end do
808 end do
809 end do
810 end do
811 call restart_close(restart, coulomb_int_file)
812
813 if (debug%info) then
814 message(1) = "Debug: Writing Coulomb integrals restart done."
815 call messages_info(1)
816 end if
817
819 end subroutine lda_u_dump_coulomb_integrals
820
821end module lda_u_io_oct_m
character(len=1), dimension(0:3), parameter, public l_notation
type(debug_t), save, public debug
Definition: debug.F90:154
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_zero
Definition: global.F90:187
complex(real64), parameter, public m_z0
Definition: global.F90:197
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:468
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:354
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:395
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:352
subroutine, public lda_u_dump(restart, namespace, this, st, mesh, ierr)
Definition: lda_u_io.F90:646
subroutine, public lda_u_write_u(this, iunit, namespace)
Definition: lda_u_io.F90:530
subroutine, public lda_u_load(restart, this, st, dftu_energy, ierr, occ_only, u_only)
Definition: lda_u_io.F90:731
subroutine lda_u_dump_coulomb_integrals(this, namespace, restart, st, mesh, ierr)
Definition: lda_u_io.F90:835
subroutine, public lda_u_write_magnetization(dir, this, ions, mesh, st, namespace)
Definition: lda_u_io.F90:483
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:579
subroutine, public lda_u_get_effectiveu(this, Ueff)
Definition: lda_u.F90:918
subroutine, public lda_u_set_effectiveu(this, Ueff)
Definition: lda_u.F90:902
subroutine, public dlda_u_get_occupations(this, occ)
Definition: lda_u.F90:3502
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:4222
integer, parameter, public dft_u_empirical
Definition: lda_u.F90:200
subroutine, public compute_acbno_u_kanamori(this, st, kanamori)
Definition: lda_u.F90:865
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:4296
subroutine, public dlda_u_set_occupations(this, occ)
Definition: lda_u.F90:3445
subroutine, public lda_u_get_effectivev(this, Veff)
Definition: lda_u.F90:952
subroutine, public lda_u_set_effectivev(this, Veff)
Definition: lda_u.F90:934
integer, parameter, public dft_u_acbn0
Definition: lda_u.F90:200
subroutine, public zlda_u_get_occupations(this, occ)
Definition: lda_u.F90:5540
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:2235
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:2309
subroutine, public zlda_u_set_occupations(this, occ)
Definition: lda_u.F90:5483
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
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:950
subroutine, public restart_write(restart, iunit, lines, nlines, ierr)
Definition: restart.F90:906
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:967
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:859
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:213
Describes mesh distribution to nodes.
Definition: mesh.F90:186
The states_elec_t class contains all electronic wave functions.
int true(void)