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