Octopus
states_elec_restart.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2016 M. Marques, A. Castro, A. Rubio, G. Bertsch, X. Andrade
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22
25 use debug_oct_m
27 use global_oct_m
29 use, intrinsic :: iso_fortran_env
33 use loct_oct_m
34 use math_oct_m
35 use mesh_oct_m
38 use mpi_oct_m
42 use parser_oct_m
45 use smear_oct_m
46 use space_oct_m
51 use string_oct_m
52 use types_oct_m
53
54 implicit none
55
56
57 private
58 public :: &
70
71contains
72
73 ! ---------------------------------------------------------
74 subroutine states_elec_look_and_load(restart, namespace, space, st, mesh, kpoints, fixed_occ, is_complex, packed)
75 type(restart_t), intent(in) :: restart
76 type(namespace_t), intent(in) :: namespace
77 class(space_t), intent(in) :: space
78 type(states_elec_t), target, intent(inout) :: st
79 class(mesh_t), intent(in) :: mesh
80 type(kpoints_t), intent(in) :: kpoints
81 logical, intent(in) :: fixed_occ
82 logical, optional, intent(in) :: is_complex
83 logical, optional, intent(in) :: packed
84
85 integer :: nkpt, dim, nst, ierr
86 real(real64), allocatable :: new_occ(:,:)
87
89
90 !check how many wfs we have
91 call states_elec_look(restart, nkpt, dim, nst, ierr)
92 if (ierr /= 0) then
93 message(1) = "Unable to read states information."
94 call messages_fatal(1, namespace=namespace)
95 end if
96
97 if (st%parallel_in_states) then
98 message(1) = "Internal error: cannot use states_elec_look_and_load when parallel in states."
99 call messages_fatal(1, namespace=namespace)
100 end if
101
102 ! Resize st%occ, retaining information there
103 allocate(new_occ(1:nst, 1:st%nik))
104 new_occ(:,:) = m_zero
105 new_occ(1:min(nst, st%nst),:) = st%occ(1:min(nst, st%nst),:)
106 safe_deallocate_a(st%occ)
107 call move_alloc(new_occ, st%occ)
108
109 ! FIXME: This wrong, one cannot just change the number of states
110 ! without updating the internal structures, in the case of parallel in states.
111 ! I think it is right now without state parallelism.
112 st%nst = nst
113 st%st_start = 1
114 st%st_end = nst
115 st%lnst = nst
117 safe_deallocate_a(st%node)
118 safe_allocate(st%node(1:st%nst))
119 st%node(:) = 0
120
121 safe_deallocate_a(st%eigenval)
122 safe_allocate(st%eigenval(1:st%nst, 1:st%nik))
123 st%eigenval = huge(st%eigenval)
124
125 if (present(is_complex)) then
126 if (is_complex) then
127 call states_elec_allocate_wfns(st, mesh, type_cmplx, packed=optional_default(packed, .false.))
128 else
129 call states_elec_allocate_wfns(st, mesh, type_float, packed=optional_default(packed, .false.))
130 end if
131 else
132 ! allow states_elec_allocate_wfns to decide for itself whether complex or real needed
133 call states_elec_allocate_wfns(st, mesh, packed=optional_default(packed, .false.))
134 end if
135
136 if (st%d%ispin == spinors) then
137 safe_allocate(st%spin(1:3, 1:st%nst, 1:st%nik))
138 st%spin = m_zero
139 end if
140
141 ! load wavefunctions
142 call states_elec_load(restart, namespace, space, st, mesh, kpoints, fixed_occ, ierr)
143 if (ierr /= 0) then
144 message(1) = "Unable to read wavefunctions."
145 call messages_fatal(1, namespace=namespace)
146 end if
147
149 end subroutine states_elec_look_and_load
150
151
152 ! ---------------------------------------------------------
153 subroutine states_elec_dump(restart, space, st, mesh, kpoints, ierr, iter, lr, st_start_writing, verbose)
154 type(restart_t), intent(in) :: restart
155 class(space_t), intent(in) :: space
156 type(states_elec_t), intent(in) :: st
157 class(mesh_t), intent(in) :: mesh
158 type(kpoints_t), intent(in) :: kpoints
159 integer, intent(out) :: ierr
160 integer, optional, intent(in) :: iter
162 type(lr_t), optional, intent(in) :: lr
163 integer, optional, intent(in) :: st_start_writing
164 logical, optional, intent(in) :: verbose
165
166 integer :: iunit_wfns, iunit_occs, iunit_states
167 integer :: err, err2(2), ik, idir, ist, idim, itot
168 integer :: root(1:p_strategy_max)
169 character(len=MAX_PATH_LEN) :: filename
170 character(len=300) :: lines(3)
171 logical :: lr_wfns_are_associated, should_write, verbose_
172 real(real64) :: kpoint(space%dim)
173 real(real64), allocatable :: dpsi(:)
174 complex(real64), allocatable :: zpsi(:)
175
176 push_sub(states_elec_dump)
177
178 verbose_ = optional_default(verbose, .true.)
179
180 ierr = 0
181
182 if (restart_skip(restart)) then
183 pop_sub(states_elec_dump)
184 return
185 end if
186
187 if (verbose_) then
188 message(1) = "Info: Writing states."
189 call print_date(trim(message(1))//' ')
190 end if
191
192 call profiling_in("RESTART_WRITE")
193
194 if (present(lr)) then
195 lr_wfns_are_associated = (allocated(lr%ddl_psi) .and. states_are_real(st)) .or. &
196 (allocated(lr%zdl_psi) .and. states_are_complex(st))
197 assert(lr_wfns_are_associated)
198 end if
199
201
202
203 iunit_states = restart_open(restart, 'states')
204 write(lines(1), '(a20,1i10)') 'nst= ', st%nst
205 write(lines(2), '(a20,1i10)') 'dim= ', st%d%dim
206 write(lines(3), '(a20,1i10)') 'nik= ', st%nik
207 call restart_write(restart, iunit_states, lines, 3, err)
208 if (err /= 0) ierr = ierr + 1
209 call restart_close(restart, iunit_states)
210
211
212 iunit_wfns = restart_open(restart, 'wfns')
213 lines(1) = '# #k-point #st #dim filename'
214 if (states_are_real(st)) then
215 lines(2) = '%Real_Wavefunctions'
216 else
217 lines(2) = '%Complex_Wavefunctions'
218 end if
219 call restart_write(restart, iunit_wfns, lines, 2, err)
220 if (err /= 0) ierr = ierr + 2
221
222
223 iunit_occs = restart_open(restart, 'occs')
224 lines(1) = '# occupations | eigenvalue[a.u.] | Im(eigenvalue) [a.u.] | k-points | k-weights | filename | ik | ist | idim'
225 lines(2) = '%Occupations_Eigenvalues_K-Points'
226 call restart_write(restart, iunit_occs, lines, 2, err)
227 if (err /= 0) ierr = ierr + 4
228
229
230 if (states_are_real(st)) then
231 safe_allocate(dpsi(1:mesh%np))
232 else
233 safe_allocate(zpsi(1:mesh%np))
234 end if
235
236 itot = 1
237 root = -1
238 err2 = 0
239 do ik = 1, st%nik
240 kpoint(1:space%dim) = &
241 kpoints%get_point(st%d%get_kpoint_index(ik), absolute_coordinates = .true.)
242
243 do ist = 1, st%nst
244 do idim = 1, st%d%dim
245 root(p_strategy_domains) = mod(itot - 1, mesh%mpi_grp%size)
246 write(filename,'(i10.10)') itot
247
248 write(lines(1), '(i8,a,i8,a,i8,3a)') ik, ' | ', ist, ' | ', idim, ' | "', trim(filename), '"'
249 call restart_write(restart, iunit_wfns, lines, 1, err)
250 if (err /= 0) err2(1) = err2(1) + 1
251
252 write(lines(1), '(e23.16,a,e23.16)') st%occ(ist,ik), ' | ', st%eigenval(ist, ik)
253 write(lines(1), '(a,a,e23.16)') trim(lines(1)), ' | ', m_zero
254 do idir = 1, space%dim
255 write(lines(1), '(a,a,e23.16)') trim(lines(1)), ' | ', kpoint(idir)
256 end do
257 write(lines(1), '(a,a,e23.16,a,i10.10,3(a,i8))') trim(lines(1)), &
258 ' | ', st%kweights(ik), ' | ', itot, ' | ', ik, ' | ', ist, ' | ', idim
259 call restart_write(restart, iunit_occs, lines, 1, err)
260 if (err /= 0) err2(1) = err2(1) + 1
261
262 should_write = st%st_start <= ist .and. ist <= st%st_end
263 if (should_write .and. present(st_start_writing)) then
264 if (ist < st_start_writing) should_write = .false.
265 end if
266
267 if (should_write) then
268 if (.not. present(lr)) then
269 if (st%d%kpt%start <= ik .and. ik <= st%d%kpt%end) then
270 if (states_are_real(st)) then
271 call states_elec_get_state(st, mesh, idim, ist, ik, dpsi)
272 call restart_write_mesh_function(restart, space, filename, mesh, dpsi, err, root = root)
273 else
274 call states_elec_get_state(st, mesh, idim, ist, ik, zpsi)
275 call restart_write_mesh_function(restart, space, filename, mesh, zpsi, err, root = root)
276 end if
277 else
278 err = 0
279 end if
280 else
281 if (st%d%kpt%start <= ik .and. ik <= st%d%kpt%end) then
282 if (states_are_real(st)) then
283 call restart_write_mesh_function(restart, space, filename, mesh, &
284 lr%ddl_psi(:, idim, ist, ik), err, root = root)
285 else
286 call restart_write_mesh_function(restart, space, filename, mesh, &
287 lr%zdl_psi(:, idim, ist, ik), err, root = root)
288 end if
289 else
290 err = 0
291 end if
292 end if
293 if (err /= 0) err2(2) = err2(2) + 1
294 end if
295
296 itot = itot + 1
297 end do ! st%d%dim
298 end do ! st%nst
299 end do ! st%nik
300 if (err2(1) /= 0) ierr = ierr + 8
301 if (err2(2) /= 0) ierr = ierr + 16
302
303 safe_deallocate_a(dpsi)
304 safe_deallocate_a(zpsi)
305
306 lines(1) = '%'
307 call restart_write(restart, iunit_occs, lines, 1, err)
308 if (err /= 0) ierr = ierr + 32
309 call restart_write(restart, iunit_wfns, lines, 1, err)
310 if (err /= 0) ierr = ierr + 64
311 if (present(iter)) then
312 write(lines(1),'(a,i7)') 'Iter = ', iter
313 call restart_write(restart, iunit_wfns, lines, 1, err)
314 if (err /= 0) ierr = ierr + 128
315 end if
316
317 call restart_close(restart, iunit_wfns)
318 call restart_close(restart, iunit_occs)
319
320 if (verbose_) then
321 message(1) = "Info: Finished writing states."
322 call print_date(trim(message(1))//' ')
323 end if
324
326
327 call profiling_out("RESTART_WRITE")
328 pop_sub(states_elec_dump)
329 return
330 end subroutine states_elec_dump
331
332
333 ! ---------------------------------------------------------
338 subroutine states_elec_load(restart, namespace, space, st, mesh, kpoints, fixed_occ, &
339 ierr, iter, lr, lowest_missing, label, verbose, skip)
340 type(restart_t), intent(in) :: restart
341 type(namespace_t), intent(in) :: namespace
342 class(space_t), intent(in) :: space
343 type(states_elec_t), intent(inout) :: st
344 class(mesh_t), intent(in) :: mesh
345 type(kpoints_t), intent(in) :: kpoints
346 logical, intent(in) :: fixed_occ
347 integer, intent(out) :: ierr
348 integer, optional, intent(out) :: iter
349 type(lr_t), optional, intent(inout) :: lr
350 integer, optional, intent(out) :: lowest_missing(:, :)
351 character(len=*), optional, intent(in) :: label
352 logical, optional, intent(in) :: verbose
353 logical, optional, intent(in) :: skip(:)
354
355 integer :: states_elec_file, wfns_file, occ_file, err, ik, ist, idir, idim
356 integer :: idone, iread, ntodo
357 character(len=12) :: filename
358 character(len=1) :: char
359 logical, allocatable :: filled(:, :, :)
360 character(len=256) :: lines(3), label_
361 character(len=50) :: str
362
363 real(real64) :: my_occ, imev, my_kweight
364 logical :: read_occ, lr_allocated, verbose_
365 logical :: integral_occs
366 real(real64), allocatable :: dpsi(:)
367 complex(real64), allocatable :: zpsi(:)
368 character(len=256), allocatable :: restart_file(:, :, :)
369 logical, allocatable :: restart_file_present(:, :, :)
370 real(real64) :: kpoint(space%dim), read_kpoint(space%dim)
371
372 integer :: iread_tmp
373 integer, allocatable :: lowest_missing_tmp(:, :)
374
375 push_sub(states_elec_load)
376
377 ierr = 0
378
379 ! make sure these intent(out)`s are initialized no matter what
380 if (present(lowest_missing)) lowest_missing = 1
381 if (present(iter)) iter = 0
382
383 if (present(skip)) then
384 assert(ubound(skip, dim = 1) == st%nst)
385 end if
386
387 if (restart_skip(restart)) then
388 ierr = -1
389 pop_sub(states_elec_load)
390 return
391 end if
392
393 call profiling_in("RESTART_READ")
394
395 verbose_ = optional_default(verbose, .true.)
396
397 if (present(label)) then
398 label_ = trim(label)
399 else
400 if (present(lr)) then
401 label_ = " for linear response"
402 else
403 label_ = ""
404 end if
405 end if
406
407 message(1) = 'Info: Reading states'
408 if (len(trim(label_)) > 0) then
409 message(1) = trim(message(1)) // trim(label_)
410 end if
411 message(1) = trim(message(1)) // "."
412 if (verbose_) call print_date(trim(message(1))//' ')
413
414 if (.not. present(lr)) then
415 st%fromScratch = .false. ! obviously, we are using restart info
416 end if
417
418 ! If one restarts a GS calculation changing the %Occupations block, one
419 ! cannot read the occupations, otherwise these overwrite the ones from
420 ! the input file. restart_fixed_occ makes that we do use the ones in the file.
421 integral_occs = .true. ! only used if restart_fixed_occ
422 if (fixed_occ) then
423 read_occ = .true.
424 st%fixed_occ = .true.
425 else
426 read_occ = .not. st%fixed_occ
427 end if
428
429 if (.not. present(lr)) then
430 st%eigenval(:, :) = m_zero
431 ! to be filled in from reading afterward
432 end if
433
434 if (.not. present(lr) .and. read_occ) then
435 st%occ(:, :) = m_zero
436 ! to be filled in from reading afterward
437 end if
438
439 ! sanity check
440 if (present(lr)) then
441 lr_allocated = (allocated(lr%ddl_psi) .and. states_are_real(st)) .or. &
442 (allocated(lr%zdl_psi) .and. states_are_complex(st))
443 assert(lr_allocated)
444 end if
445
446 states_elec_file = restart_open(restart, 'states')
447 ! sanity check on spin/k-points. Example file 'states':
448 ! nst= 2
449 ! dim= 1
450 ! nik= 2
451 call restart_read(restart, states_elec_file, lines, 3, err)
452 if (err /= 0) then
453 ierr = ierr - 2
454 else
455 read(lines(2), *) str, idim
456 read(lines(3), *) str, ik
457 if (idim == 2 .and. st%d%dim == 1) then
458 write(message(1),'(a)') 'Incompatible restart information: saved calculation is spinors, this one is not.'
459 call messages_warning(1, namespace=namespace)
460 ierr = ierr - 2**2
461 end if
462 if (idim == 1 .and. st%d%dim == 2) then
463 write(message(1),'(a)') 'Incompatible restart information: this calculation is spinors, saved one is not.'
464 call messages_warning(1, namespace=namespace)
465 ierr = ierr - 2**3
466 end if
467 if (ik < st%nik) then
468 write(message(1),'(a)') 'Incompatible restart information: not enough k-points.'
469 write(message(2),'(2(a,i6))') 'Expected ', st%nik, ' > Read ', ik
470 call messages_warning(2, namespace=namespace)
471 end if
472 ! We will check that they are the right k-points later, so we do not need to put a specific error here.
473 end if
474 call restart_close(restart, states_elec_file)
475
476
477 ! open files to read
478 wfns_file = restart_open(restart, 'wfns')
479 occ_file = restart_open(restart, 'occs')
480 call restart_read(restart, wfns_file, lines, 2, err)
481 if (err /= 0) then
482 ierr = ierr - 2**5
483 else if (states_are_real(st)) then
484 read(lines(2), '(a)') str
485 if (str(2:8) == 'Complex') then
486 message(1) = "Cannot read real states from complex wavefunctions."
487 call messages_warning(1, namespace=namespace)
488 ierr = ierr - 2**6
489 else if (str(2:5) /= 'Real') then
490 message(1) = "Restart file 'wfns' does not specify real/complex; cannot check compatibility."
491 call messages_warning(1, namespace=namespace)
492 end if
493 end if
494 ! complex can be restarted from real, so there is no problem.
495
496 ! Skip two lines.
497 call restart_read(restart, occ_file, lines, 2, err)
498 if (err /= 0) ierr = ierr - 2**7
499
500 ! If any error occured up to this point then it is not worth continuing,
501 ! as there something fundamentally wrong with the restart files
502 if (ierr /= 0) then
503 call restart_close(restart, wfns_file)
504 call restart_close(restart, occ_file)
505 call profiling_out("RESTART_READ")
506 pop_sub(states_elec_load)
507 return
508 end if
509
510 if (states_are_real(st)) then
511 safe_allocate(dpsi(1:mesh%np))
512 else
513 safe_allocate(zpsi(1:mesh%np))
514 end if
515
516 safe_allocate(restart_file(1:st%d%dim, st%st_start:st%st_end, 1:st%nik))
517 safe_allocate(restart_file_present(1:st%d%dim, st%st_start:st%st_end, 1:st%nik))
518 restart_file_present = .false.
519
520 ! Next we read the list of states from the files.
521 ! Errors in reading the information of a specific state from the files are ignored
522 ! at this point, because later we will skip reading the wavefunction of that state.
523 do
524 call restart_read(restart, wfns_file, lines, 1, err)
525 if (err == 0) then
526 read(lines(1), '(a)') char
527 if (char == '%') then
528 !We reached the end of the file
529 exit
530 else
531 read(lines(1), *) ik, char, ist, char, idim, char, filename
532 end if
533 end if
534
535 if (err /= 0 .or. index_is_wrong()) then
536 call restart_read(restart, occ_file, lines, 1, err)
537 cycle
538 end if
539
540 if (ist >= st%st_start .and. ist <= st%st_end .and. &
541 ik >= st%d%kpt%start .and. ik <= st%d%kpt%end) then
542
543 restart_file(idim, ist, ik) = trim(filename)
544 restart_file_present(idim, ist, ik) = .true.
545 end if
546
547 call restart_read(restart, occ_file, lines, 1, err)
548 if (.not. present(lr)) then ! do not read eigenvalues or occupations when reading linear response
549 ! # occupations | eigenvalue[a.u.] | Im(eigenvalue) [a.u.] | k-points | k-weights | filename | ik | ist | idim
550
551 if (err == 0) then
552 read(lines(1), *) my_occ, char, st%eigenval(ist, ik), char, imev, char, &
553 (read_kpoint(idir), char, idir = 1, space%dim), my_kweight
554 ! we do not want to read the k-weights, we have already set them appropriately
555 else
556 ! There is a problem with this states information, so we skip it.
557 if (ist >= st%st_start .and. ist <= st%st_end .and. &
558 ik >= st%d%kpt%start .and. ik <= st%d%kpt%end) then
559 restart_file_present(idim, ist, ik) = .false.
560 end if
561 cycle
562 end if
563
564 kpoint(1:space%dim) = &
565 kpoints%get_point(st%d%get_kpoint_index(ik), absolute_coordinates = .true.)
566 ! FIXME: maybe should ignore ik and just try to match actual vector k-points?
567 if (any(abs(kpoint(1:space%dim) - read_kpoint(1:space%dim)) > 1e-12_real64)) then
568 ! write only once for each k-point so as not to be too verbose
569 if (ist == 1) then
570 write(message(1),'(a,i6)') 'Incompatible restart information: k-point mismatch for ik ', ik
571 write(message(2),'(a,99f18.12)') ' Expected : ', kpoint(1:space%dim)
572 write(message(3),'(a,99f18.12)') ' Read : ', read_kpoint(1:space%dim)
573 call messages_warning(3, namespace=namespace)
574 end if
575 if (ist >= st%st_start .and. ist <= st%st_end .and. &
576 ik >= st%d%kpt%start .and. ik <= st%d%kpt%end) then
577 restart_file_present(idim, ist, ik) = .false.
578 end if
579 end if
580
581 if (read_occ) then
582 st%occ(ist, ik) = my_occ
583 integral_occs = integral_occs .and. &
584 abs((st%occ(ist, ik) - st%smear%el_per_state) * st%occ(ist, ik)) <= m_epsilon
585 end if
586 end if
587 end do
588
589 if (present(iter)) then
590 call restart_read(restart, wfns_file, lines, 1, err)
591 if (err /= 0) then
592 ierr = ierr - 2**8
593 else
594 read(lines(1), *) filename, filename, iter
595 end if
596 end if
597
598 call restart_close(restart, wfns_file)
599 call restart_close(restart, occ_file)
600
601 ! Now we read the wavefunctions. At this point we need to have all the information from the
602 ! states, occs, and wfns files in order to avoid serialisation of reads, as restart_read
603 ! works as a barrier.
604
605 safe_allocate(filled(1:st%d%dim, st%st_start:st%st_end, st%d%kpt%start:st%d%kpt%end))
606 filled = .false.
607
608 if (present(lowest_missing)) lowest_missing = st%nst + 1
609
610 iread = 0
611 if (mpi_grp_is_root(mpi_world) .and. verbose_) then
612 idone = 0
613 ntodo = st%lnst*st%d%kpt%nlocal*st%d%dim
614 call loct_progress_bar(-1, ntodo)
615 end if
616
617 do ik = st%d%kpt%start, st%d%kpt%end
618 do ist = st%st_start, st%st_end
619 if (present(skip)) then
620 if (skip(ist)) cycle
621 end if
622
623 do idim = 1, st%d%dim
624
625 if (.not. restart_file_present(idim, ist, ik)) then
626 if (present(lowest_missing)) then
627 lowest_missing(idim, ik) = min(lowest_missing(idim, ik), ist)
628 end if
629 cycle
630 end if
631
632 if (states_are_real(st)) then
633 call restart_read_mesh_function(restart, space, restart_file(idim, ist, ik), mesh, dpsi, err)
634 else
635 call restart_read_mesh_function(restart, space, restart_file(idim, ist, ik), mesh, zpsi, err)
636 end if
637
638 if (states_are_real(st)) then
639 if (.not. present(lr)) then
640 call states_elec_set_state(st, mesh, idim, ist, ik, dpsi)
641 else
642 call lalg_copy(mesh%np, dpsi, lr%ddl_psi(:, idim, ist, ik))
643 end if
644 else
645 if (.not. present(lr)) then
646 call states_elec_set_state(st, mesh, idim, ist, ik, zpsi)
647 else
648 call lalg_copy(mesh%np, zpsi, lr%zdl_psi(:, idim, ist, ik))
649 end if
650 end if
651
652
653 if (err == 0) then
654 filled(idim, ist, ik) = .true.
655 iread = iread + 1
656 else if (present(lowest_missing)) then
657 lowest_missing(idim, ik) = min(lowest_missing(idim, ik), ist)
658 end if
659
660 if (mpi_grp_is_root(mpi_world) .and. verbose_) then
661 idone = idone + 1
662 call loct_progress_bar(idone, ntodo)
663 end if
664
665 end do
666 end do
667 end do
668
669 safe_deallocate_a(dpsi)
670 safe_deallocate_a(zpsi)
671 safe_deallocate_a(restart_file)
672 safe_deallocate_a(restart_file_present)
673
674 if (mpi_grp_is_root(mpi_world) .and. verbose_) then
675 call messages_new_line()
676 end if
677
678 if (st%parallel_in_states .or. st%d%kpt%parallel) then
679 iread_tmp = iread
680 call st%st_kpt_mpi_grp%allreduce(iread_tmp, iread, 1, mpi_integer, mpi_sum)
681 end if
682
683 if (st%d%kpt%parallel) then
684 ! make sure all tasks know lowest_missing over all k-points
685 if (present(lowest_missing)) then
686 safe_allocate(lowest_missing_tmp(1:st%d%dim, 1:st%nik))
687 lowest_missing_tmp = lowest_missing
688 call st%d%kpt%mpi_grp%allreduce(lowest_missing_tmp(1,1), lowest_missing(1,1), st%d%dim*st%nik, &
689 mpi_integer, mpi_min)
690 safe_deallocate_a(lowest_missing_tmp)
691 end if
692 end if
693
694 if (fixed_occ .and. iread == st%nst * st%nik * st%d%dim) then
695 ! reset to overwrite whatever smearing may have been set earlier
696 ! only do this if all states could be loaded from restart files
697 call smear_init(st%smear, namespace, st%d%ispin, fixed_occ = .true., integral_occs = integral_occs, kpoints = kpoints)
698 end if
699
700 if (.not. present(lr) .and. .not. present(skip)) call fill_random()
701 ! it is better to initialize lr wfns to zero
702
703 safe_deallocate_a(filled)
704
705 if (ierr == 0 .and. iread /= st%nst * st%nik * st%d%dim) then
706 if (iread > 0) then
707 ierr = iread
708 else
709 ierr = -1
710 end if
711 ! otherwise ierr = 0 would mean either all was read correctly, or nothing at all was read!
712
713 if (.not. present(lr)) then
714 write(str, '(a,i5)') 'Reading states.'
715 else
716 write(str, '(a,i5)') 'Reading states information for linear response.'
717 end if
718 call messages_print_with_emphasis(msg=trim(str), namespace=namespace)
719 if (.not. present(skip)) then
720 write(message(1),'(a,i6,a,i6,a)') 'Only ', iread,' files out of ', &
721 st%nst * st%nik * st%d%dim, ' could be read.'
722 else
723 write(message(1),'(a,i6,a,i6,a)') 'Only ', iread,' files out of ', &
724 st%nst * st%nik * st%d%dim, ' were loaded.'
725 ierr = 0
726 end if
727 call messages_info(1, namespace=namespace)
728 call messages_print_with_emphasis(namespace=namespace)
729 end if
730
731 message(1) = 'Info: States reading done.'
732 if (verbose_) call print_date(trim(message(1))//' ')
733
734 call profiling_out("RESTART_READ")
735 pop_sub(states_elec_load)
736
737 contains
738
739 ! ---------------------------------------------------------
740 subroutine fill_random() !< put random function in orbitals that could not be read.
742
743 do ik = st%d%kpt%start, st%d%kpt%end
744
745 do ist = st%st_start, st%st_end
746 do idim = 1, st%d%dim
747 if (filled(idim, ist, ik)) cycle
748
749 call states_elec_generate_random(st, mesh, kpoints, ist, ist, ik, ik)
750 end do
751 end do
752 end do
753
755 end subroutine fill_random
756
757 ! ---------------------------------------------------------
758
759 logical function index_is_wrong() !< .true. if the index (idim, ist, ik) is not present in st structure...
761
762 if (idim > st%d%dim .or. idim < 1 .or. &
763 ist > st%nst .or. ist < 1 .or. &
764 ik > st%nik .or. ik < 1) then
766 else
767 index_is_wrong = .false.
768 end if
769
771 end function index_is_wrong
772
773 end subroutine states_elec_load
774
775
776 subroutine states_elec_dump_rho(restart, space, st, mesh, ierr, iter)
777 type(restart_t), intent(in) :: restart
778 class(space_t), intent(in) :: space
779 type(states_elec_t), intent(in) :: st
780 class(mesh_t), intent(in) :: mesh
781 integer, intent(out) :: ierr
782 integer, optional, intent(in) :: iter
783
784 integer :: iunit, isp, err, err2(2)
785 character(len=MAX_PATH_LEN) :: filename
786 character(len=300) :: lines(2)
787
788 push_sub(states_elec_dump_rho)
789
790 ierr = 0
791
792 if (restart_skip(restart)) then
793 pop_sub(states_elec_dump_rho)
794 return
795 end if
796
797 message(1) = "Debug: Writing density restart."
798 call messages_info(1, debug_only=.true.)
799
801
802 !write the densities
803 iunit = restart_open(restart, 'density')
804 lines(1) = '# #spin #nspin filename'
805 lines(2) = '%densities'
806 call restart_write(restart, iunit, lines, 2, err)
807 if (err /= 0) ierr = ierr + 1
808
809 err2 = 0
810 do isp = 1, st%d%nspin
811 filename = get_filename_with_spin('density', st%d%nspin, isp)
812 write(lines(1), '(i8,a,i8,a)') isp, ' | ', st%d%nspin, ' | "'//trim(adjustl(filename))//'"'
813 call restart_write(restart, iunit, lines, 1, err)
814 if (err /= 0) err2(1) = err2(1) + 1
815
816 call restart_write_mesh_function(restart, space, filename, mesh, st%rho(:,isp), err)
817 if (err /= 0) err2(2) = err2(2) + 1
818
819 end do
820 if (err2(1) /= 0) ierr = ierr + 2
821 if (err2(2) /= 0) ierr = ierr + 4
822
823 lines(1) = '%'
824 call restart_write(restart, iunit, lines, 1, err)
825 if (err /= 0) ierr = ierr + 8
826 if (present(iter)) then
827 write(lines(1),'(a,i7)') 'Iter = ', iter
828 call restart_write(restart, iunit, lines, 1, err)
829 if (err /= 0) ierr = ierr + 16
830 end if
831 call restart_close(restart, iunit)
832
834
835 message(1) = "Debug: Writing density restart done."
836 call messages_info(1, debug_only=.true.)
837
838 pop_sub(states_elec_dump_rho)
839 end subroutine states_elec_dump_rho
840
841
842 ! ---------------------------------------------------------
843 subroutine states_elec_load_rho(restart, space, st, mesh, ierr)
844 type(restart_t), intent(in) :: restart
845 class(space_t), intent(in) :: space
846 type(states_elec_t), intent(inout) :: st
847 class(mesh_t), intent(in) :: mesh
848 integer, intent(out) :: ierr
849
850 integer :: err, err2, isp
851 character(len=MAX_PATH_LEN) :: filename
853 push_sub(states_elec_load_rho)
854
855 ierr = 0
856
857 if (restart_skip(restart)) then
858 ierr = -1
859 pop_sub(states_elec_load_rho)
860 return
861 end if
862
863 message(1) = "Debug: Reading density restart."
864 call messages_info(1, debug_only=.true.)
865
866 ! skip for now, since we know what the files are going to be called
867 !read the densities
868! iunit_rho = io_open(trim(dir)//'/density', restart%namespace, action='write')
869! call iopar_read(st%dom_st_kpt_mpi_grp, iunit_rho, line, err)
870! call iopar_read(st%dom_st_kpt_mpi_grp, iunit_rho, line, err)
871! we could read the iteration 'iter' too, not sure if that is useful.
872
873 err2 = 0
874 do isp = 1, st%d%nspin
875 filename = get_filename_with_spin('density', st%d%nspin, isp)
876! if (mpi_grp_is_root(st%dom_st_kpt_mpi_grp)) then
877! read(iunit_rho, '(i8,a,i8,a)') isp, ' | ', st%d%nspin, ' | "'//trim(adjustl(filename))//'"'
878! end if
879 call restart_read_mesh_function(restart, space, filename, mesh, st%rho(:,isp), err)
880 if (err /= 0) err2 = err2 + 1
881
882 end do
883 if (err2 /= 0) ierr = ierr + 1
884
885 message(1) = "Debug: Reading density restart done."
886 call messages_info(1, debug_only=.true.)
887
888 pop_sub(states_elec_load_rho)
889 end subroutine states_elec_load_rho
890
891 subroutine states_elec_dump_frozen(restart, space, st, mesh, ierr)
892 type(restart_t), intent(in) :: restart
893 class(space_t), intent(in) :: space
894 type(states_elec_t), intent(in) :: st
895 class(mesh_t), intent(in) :: mesh
896 integer, intent(out) :: ierr
897
898 integer :: isp, err, err2(2), idir
899 character(len=MAX_PATH_LEN) :: filename
900
902
903 ierr = 0
904
905 assert(allocated(st%frozen_rho))
906
907 if (restart_skip(restart)) then
909 return
910 end if
911
912 message(1) = "Debug: Writing frozen densities restart."
913 call messages_info(1, debug_only=.true.)
914
916
917 err2 = 0
918 do isp = 1, st%d%nspin
919 filename = get_filename_with_spin('frozen_rho', st%d%nspin, isp)
920
921 call restart_write_mesh_function(restart, space, filename, mesh, st%frozen_rho(:,isp), err)
922 if (err /= 0) err2(2) = err2(2) + 1
923
924 if (allocated(st%frozen_tau)) then
925 filename = get_filename_with_spin('frozen_tau', st%d%nspin, isp)
926 call restart_write_mesh_function(restart, space, filename, mesh, st%frozen_tau(:,isp), err)
927 if (err /= 0) err2 = err2 + 1
928 end if
929
930 if (allocated(st%frozen_gdens)) then
931 do idir = 1, space%dim
932 if (st%d%nspin == 1) then
933 write(filename, fmt='(a,i1)') 'frozen_gdens-dir', idir
934 else
935 write(filename, fmt='(a,i1,a,i1)') 'frozen_gdens-dir', idir, '-', isp
936 end if
937 call restart_write_mesh_function(restart, space, filename, mesh, st%frozen_gdens(:,idir,isp), err)
938 if (err /= 0) err2 = err2 + 1
939 end do
940 end if
941
942 if (allocated(st%frozen_ldens)) then
943 filename = get_filename_with_spin('frozen_ldens', st%d%nspin, isp)
944 call restart_write_mesh_function(restart, space, filename, mesh, st%frozen_ldens(:,isp), err)
945 if (err /= 0) err2 = err2 + 1
946 end if
947
948 end do
949 if (err2(1) /= 0) ierr = ierr + 2
950 if (err2(2) /= 0) ierr = ierr + 4
951
953
954 message(1) = "Debug: Writing frozen densities restart done."
955 call messages_info(1, debug_only=.true.)
956
958 end subroutine states_elec_dump_frozen
959
960
961 ! ---------------------------------------------------------
962 subroutine states_elec_load_frozen(restart, space, st, mesh, ierr)
963 type(restart_t), intent(in) :: restart
964 class(space_t), intent(in) :: space
965 type(states_elec_t), intent(inout) :: st
966 class(mesh_t), intent(in) :: mesh
967 integer, intent(out) :: ierr
968
969 integer :: err, err2, isp, idir
970 character(len=MAX_PATH_LEN) :: filename
971
973
974 assert(allocated(st%frozen_rho))
975
976 ierr = 0
977
978 if (restart_skip(restart)) then
979 ierr = -1
981 return
982 end if
983
984 message(1) = "Debug: Reading densities restart."
985 call messages_info(1, debug_only=.true.)
986
987 err2 = 0
988 do isp = 1, st%d%nspin
989 filename = get_filename_with_spin('frozen_rho', st%d%nspin, isp)
990 call restart_read_mesh_function(restart, space, filename, mesh, st%frozen_rho(:,isp), err)
991 if (err /= 0) err2 = err2 + 1
992
993 if (allocated(st%frozen_tau)) then
994 filename = get_filename_with_spin('frozen_tau', st%d%nspin, isp)
995 call restart_read_mesh_function(restart, space, filename, mesh, st%frozen_tau(:,isp), err)
996 if (err /= 0) err2 = err2 + 1
997 end if
998
999 if (allocated(st%frozen_gdens)) then
1000 do idir = 1, space%dim
1001 if (st%d%nspin == 1) then
1002 write(filename, fmt='(a,i1)') 'frozen_gdens-dir', idir
1003 else
1004 write(filename, fmt='(a,i1,a,i1)') 'frozen_gdens-dir', idir, '-', isp
1005 end if
1006 call restart_read_mesh_function(restart, space, filename, mesh, st%frozen_gdens(:,idir,isp), err)
1007 if (err /= 0) err2 = err2 + 1
1008 end do
1009 end if
1010
1011 if (allocated(st%frozen_ldens)) then
1012 filename = get_filename_with_spin('frozen_ldens', st%d%nspin, isp)
1013 call restart_read_mesh_function(restart, space, filename, mesh, st%frozen_ldens(:,isp), err)
1014 if (err /= 0) err2 = err2 + 1
1015 end if
1016
1017 end do
1018 if (err2 /= 0) ierr = ierr + 1
1019
1020 message(1) = "Debug: Reading frozen densities restart done."
1021 call messages_info(1, debug_only=.true.)
1022
1024 end subroutine states_elec_load_frozen
1025
1026
1027 ! ---------------------------------------------------------
1030 subroutine states_elec_read_user_def_orbitals(mesh, namespace, space, st)
1031 class(mesh_t), intent(in) :: mesh
1032 type(namespace_t), intent(in) :: namespace
1033 class(space_t), intent(in) :: space
1034 type(states_elec_t), intent(inout) :: st
1035
1036 type(block_t) :: blk
1037 integer :: ip, id, is, ik, nstates, state_from, ierr, ncols
1038 integer :: ib, idim, inst, inik, normalize
1039 real(real64) :: xx(space%dim), rr, psi_re, psi_im
1040 character(len=150) :: filename
1041 complex(real64), allocatable :: zpsi(:, :)
1042
1043 integer, parameter :: &
1044 state_from_formula = 1, &
1045 state_from_file = -10010, &
1046 normalize_yes = 1, &
1047 normalize_no = 0
1048
1050
1051 !%Variable UserDefinedStates
1052 !%Type block
1053 !%Section States
1054 !%Description
1055 !% Instead of using the ground state as initial state for
1056 !% time-propagations it might be interesting in some cases
1057 !% to specify alternate states. Like with user-defined
1058 !% potentials, this block allows you to specify formulas for
1059 !% the orbitals at <i>t</i>=0.
1060 !%
1061 !% Example:
1062 !%
1063 !% <tt>%UserDefinedStates
1064 !% <br>&nbsp;&nbsp; 1 | 1 | 1 | formula | "exp(-r^2)*exp(-i*0.2*x)" | normalize_yes
1065 !% <br>%</tt>
1066 !%
1067 !% The first column specifies the component of the spinor,
1068 !% the second column the number of the state and the third
1069 !% contains <i>k</i>-point and spin quantum numbers. Column four
1070 !% indicates that column five should be interpreted as a formula
1071 !% for the corresponding orbital.
1072 !%
1073 !% Alternatively, if column four states <tt>file</tt> the state will
1074 !% be read from the file given in column five.
1075 !%
1076 !% <tt>%UserDefinedStates
1077 !% <br>&nbsp;&nbsp; 1 | 1 | 1 | file | "/path/to/file" | normalize_no
1078 !% <br>%</tt>
1079 !%
1080 !% Octopus reads first the ground-state orbitals from
1081 !% the <tt>restart/gs</tt> directory. Only the states that are
1082 !% specified in the above block will be overwritten with
1083 !% the given analytic expression for the orbital.
1084 !%
1085 !% The sixth (optional) column indicates whether <tt>Octopus</tt> should renormalize
1086 !% the orbital. The default (no sixth column given) is to renormalize.
1087 !%
1088 !%Option file -10010
1089 !% Read initial orbital from file.
1090 !% Accepted file formats, detected by extension: obf, ncdf and csv (real only).
1091 !% For csv files, the following formatting rules apply:
1092 !%
1093 !% The functions in this file read an array from an ascii matrix (csv) file.
1094 !% Format with values "valueXYZ" as follows
1095 !%
1096 !% File values.csv:
1097 !% --------
1098 !% value111 value211 value311 value411
1099 !% value121 value221 value321 value421
1100 !% value131 value231 value331 value431
1101 !%
1102 !% value112 value212 value312 value412
1103 !% value122 value222 value322 value422
1104 !% value132 value232 value332 value432
1105 !%
1106 !% value113 value213 value313 value413
1107 !% value123 value223 value323 value423
1108 !% value133 value233 value333 value433
1109 !% --------
1110 !%
1111 !% That is, every line scans the x-coordinate, every XY-plane as a table of values
1112 !% and all XY-planes separated by an empty row.
1113 !%
1114 !% The given matrix is interpolated/stretched to fit the calculation
1115 !% box defined in input file.
1116 !%
1117 !% Note that there must not be any empty line at the end of the file.
1118 !%
1119 !% Calculation box shape must be "parallelepiped".
1120 !%
1121 !% The delimiter can be a tab, a comma or a space.
1122 !%Option formula 1
1123 !% Calculate initial orbital by given analytic expression.
1124 !%Option normalize_yes 1
1125 !% Normalize orbitals (default).
1126 !%Option normalize_no 0
1127 !% Do not normalize orbitals.
1128 !%End
1129 if (parse_block(namespace, 'UserDefinedStates', blk) == 0) then
1130
1131 call messages_print_with_emphasis(msg=trim('Substitution of orbitals'), namespace=namespace)
1132
1133 ! find out how many lines (i.e. states) the block has
1134 nstates = parse_block_n(blk)
1135
1136 safe_allocate(zpsi(1:mesh%np, 1:st%d%dim))
1137
1138 ! read all lines
1139 do ib = 1, nstates
1140 ! Check that number of columns is five or six.
1141 ncols = parse_block_cols(blk, ib - 1)
1142 if (ncols < 5 .or. ncols > 6) then
1143 message(1) = 'Each line in the UserDefinedStates block must have'
1144 message(2) = 'five or six columns.'
1145 call messages_fatal(2, namespace=namespace)
1146 end if
1147
1148 call parse_block_integer(blk, ib - 1, 0, idim)
1149 call parse_block_integer(blk, ib - 1, 1, inst)
1150 call parse_block_integer(blk, ib - 1, 2, inik)
1151
1152 ! Calculate from expression or read from file?
1153 call parse_block_integer(blk, ib - 1, 3, state_from)
1154
1155 ! loop over all states
1156 do id = 1, st%d%dim
1157 do is = 1, st%nst
1158 do ik = 1, st%nik
1159
1160 if (.not.(id == idim .and. is == inst .and. ik == inik )) cycle
1161
1162 ! We first parse the value, such that the parsing and stdout looks correct
1163 ! However, if rank 0 is not responsible for the parsing, the variables in the formula
1164 ! will still appear as not parsed (pure k-point/state parallelization case).
1165 ! This could be solved by a dummy call to parse_expression
1166 select case (state_from)
1167 case (state_from_formula)
1168 ! parse formula string
1169 call parse_block_string( &
1170 blk, ib - 1, 4, st%user_def_states(id, is, ik))
1171
1172 write(message(1), '(a,3i5)') 'Substituting state of orbital with k, ist, dim = ', ik, is, id
1173 write(message(2), '(2a)') ' with the expression:'
1174 write(message(3), '(2a)') ' ',trim(st%user_def_states(id, is, ik))
1175 call messages_info(3, namespace=namespace)
1176
1177 ! convert to C string
1178 call conv_to_c_string(st%user_def_states(id, is, ik))
1179
1180 case (state_from_file)
1181 ! The input format can be coded in column four now. As it is
1182 ! not used now, we just say "file".
1183 ! Read the filename.
1184 call parse_block_string(blk, ib - 1, 4, filename)
1185
1186 write(message(1), '(a,3i5)') 'Substituting state of orbital with k, ist, dim = ', ik, is, id
1187 write(message(2), '(2a)') ' with data from file:'
1188 write(message(3), '(2a)') ' ',trim(filename)
1189 call messages_info(3, namespace=namespace)
1190 end select
1191
1192
1193 ! does the block entry match and is this node responsible?
1194 if (.not.(st%st_start <= is .and. st%st_end >= is &
1195 .and. st%d%kpt%start <= ik .and. st%d%kpt%end >= ik)) cycle
1196
1197 select case (state_from)
1198
1199 case (state_from_formula)
1200 ! fill states with user-defined formulas
1201 do ip = 1, mesh%np
1202 xx = mesh%x(ip, :)
1203 rr = norm2(xx)
1204
1205 ! parse user-defined expressions
1206 call parse_expression(psi_re, psi_im, space%dim, xx, rr, m_zero, st%user_def_states(id, is, ik))
1207 ! fill state
1208 zpsi(ip, 1) = cmplx(psi_re, psi_im, real64)
1209 end do
1210
1211 case (state_from_file)
1212 ! finally read the state
1213 call zio_function_input(filename, namespace, space, mesh, zpsi(:, 1), ierr)
1214 if (ierr > 0) then
1215 message(1) = 'Could not read the file!'
1216 write(message(2),'(a,i1)') 'Error code: ', ierr
1217 call messages_fatal(2, namespace=namespace)
1218 end if
1219
1220 case default
1221 message(1) = 'Wrong entry in UserDefinedStates, column 4.'
1222 message(2) = 'You may state "formula" or "file" here.'
1223 call messages_fatal(2, namespace=namespace)
1224 end select
1225
1226 ! normalize orbital
1227 if (parse_block_cols(blk, ib - 1) == 6) then
1228 call parse_block_integer(blk, ib - 1, 5, normalize)
1229 else
1230 normalize = 1
1231 end if
1232 select case (normalize)
1233 case (normalize_no)
1234 call states_elec_set_state(st, mesh, id, is, ik, zpsi(:, 1))
1235 case (normalize_yes)
1236 assert(st%d%dim == 1)
1237 call zmf_normalize(mesh, st%d%dim, zpsi)
1238 call states_elec_set_state(st, mesh, is, ik, zpsi)
1239 case default
1240 message(1) = 'The sixth column in UserDefinedStates may either be'
1241 message(2) = '"normalize_yes" or "normalize_no"'
1242 call messages_fatal(2, namespace=namespace)
1243 end select
1244
1245 end do
1246 end do
1247 end do
1248
1249 end do
1250
1251 safe_deallocate_a(zpsi)
1252
1253 call parse_block_end(blk)
1254 call messages_print_with_emphasis(namespace=namespace)
1255
1256 else
1257 call messages_variable_is_block(namespace, 'UserDefinedStates')
1258 end if
1259
1262
1263 ! ---------------------------------------------------------
1264 ! This is needed as for the generalized Bloch theorem we need to label
1265 ! the states from the expectation value of Sz computed from the GS.
1266 subroutine states_elec_dump_spin(restart, st, ierr)
1267 type(restart_t), intent(in) :: restart
1268 type(states_elec_t), intent(in) :: st
1269 integer, intent(out) :: ierr
1270
1271 integer :: iunit_spin
1272 integer :: err, err2(2), ik, ist
1273 character(len=300) :: lines(3)
1274
1275 push_sub(states_elec_dump_spin)
1276
1277 ierr = 0
1278
1279 if (restart_skip(restart)) then
1280 pop_sub(states_elec_dump_spin)
1281 return
1282 end if
1283
1284 call profiling_in("RESTART_WRITE_SPIN")
1285
1287
1288 iunit_spin = restart_open(restart, 'spin')
1289 lines(1) = '# #k-point #st #spin(x) spin(y) spin(z)'
1290 call restart_write(restart, iunit_spin, lines, 1, err)
1291 if (err /= 0) ierr = ierr + 1
1292
1293 err2 = 0
1294 do ik = 1, st%nik
1295 do ist = 1, st%nst
1296 write(lines(1), '(i8,a,i8,3(a,f18.12))') ik, ' | ', ist, ' | ', &
1297 st%spin(1,ist,ik), ' | ', st%spin(2,ist,ik),' | ', st%spin(3,ist,ik)
1298 call restart_write(restart, iunit_spin, lines, 1, err)
1299 if (err /= 0) err2(1) = err2(1) + 1
1300 end do ! st%nst
1301 end do ! st%nik
1302 lines(1) = '%'
1303 call restart_write(restart, iunit_spin, lines, 1, err)
1304 if (err2(1) /= 0) ierr = ierr + 8
1305 if (err2(2) /= 0) ierr = ierr + 16
1306
1307 call restart_close(restart, iunit_spin)
1308
1310
1311 call profiling_out("RESTART_WRITE_SPIN")
1312 pop_sub(states_elec_dump_spin)
1313 end subroutine states_elec_dump_spin
1314
1315
1316
1317
1318 ! ---------------------------------------------------------
1323 subroutine states_elec_load_spin(restart, st, ierr)
1324 type(restart_t), intent(in) :: restart
1325 type(states_elec_t), intent(inout) :: st
1326 integer, intent(out) :: ierr
1327
1328 integer :: spin_file, err, ik, ist
1329 character(len=256) :: lines(3)
1330 real(real64) :: spin(3)
1331 character(len=1) :: char
1332
1333
1334 push_sub(states_elec_load_spin)
1335
1336 ierr = 0
1337
1338 if (restart_skip(restart)) then
1339 ierr = -1
1340 pop_sub(states_elec_load_spin)
1341 return
1342 end if
1343
1344 call profiling_in("RESTART_READ_SPIN")
1345
1346 ! open files to read
1347 spin_file = restart_open(restart, 'spin')
1348 ! Skip two lines.
1349 call restart_read(restart, spin_file, lines, 1, err)
1350 if (err /= 0) ierr = ierr - 2**7
1351
1352 ! If any error occured up to this point then it is not worth continuing,
1353 ! as there something fundamentally wrong with the restart files
1354 if (ierr /= 0) then
1355 call restart_close(restart, spin_file)
1356 call profiling_out("RESTART_READ_SPIN")
1357 pop_sub(states_elec_load_spin)
1358 return
1359 end if
1360
1361 ! Next we read the list of states from the files.
1362 ! Errors in reading the information of a specific state from the files are ignored
1363 ! at this point, because later we will skip reading the wavefunction of that state.
1364 do
1365 call restart_read(restart, spin_file, lines, 1, err)
1366 read(lines(1), '(a)') char
1367 if (char == '%') then
1368 !We reached the end of the file
1369 exit
1370 else
1371 read(lines(1), *) ik, char, ist, char, spin(1), char, spin(2), char, spin(3)
1372 end if
1373
1374 if (err /= 0) cycle
1375
1376 st%spin(1:3, ist, ik) = spin(1:3)
1377 end do
1378
1379 call restart_close(restart, spin_file)
1380
1381 call profiling_out("RESTART_READ_SPIN")
1382 pop_sub(states_elec_load_spin)
1383 end subroutine states_elec_load_spin
1384
1385 ! ---------------------------------------------------------
1386 subroutine states_elec_transform(st, namespace, space, restart, mesh, kpoints, prefix)
1387 type(states_elec_t), intent(inout) :: st
1388 type(namespace_t), intent(in) :: namespace
1389 class(space_t), intent(in) :: space
1390 type(restart_t), intent(inout) :: restart
1391 class(mesh_t), intent(in) :: mesh
1392 type(kpoints_t), intent(in) :: kpoints
1393 character(len=*), optional, intent(in) :: prefix
1394
1395 type(states_elec_t) :: stin
1396 type(block_t) :: blk
1397 complex(real64), allocatable :: rotation_matrix(:,:), psi(:, :)
1398 integer :: ist, jst, ncols, iqn
1399 character(len=256) :: block_name
1400
1401 push_sub(states_elec_transform)
1402
1403 block_name = trim(optional_default(prefix, "")) // "TransformStates"
1404
1405 !%Variable TransformStates
1406 !%Type block
1407 !%Default no
1408 !%Section States
1409 !%Description
1410 !% Before starting the <tt>td</tt> calculation, the initial states (that are
1411 !% read from the <tt>restart/gs</tt> directory, which should have been
1412 !% generated in a previous ground-state calculation) can be "transformed"
1413 !% among themselves. The block <tt>TransformStates</tt> gives the transformation matrix
1414 !% to be used. The number of rows and columns of the matrix should equal the number
1415 !% of the states present in the time-dependent calculation (the independent
1416 !% spin and <i>k</i>-point subspaces are all transformed equally); the number of
1417 !% columns should be equal to the number of states present in the
1418 !% <tt>restart/gs</tt> directory. This number may be different: for example,
1419 !% one could have run previously in <tt>unocc</tt> mode in order to obtain unoccupied
1420 !% Kohn-Sham states, and therefore <tt>restart/gs</tt> will contain more states.
1421 !% These states can be used in the transformation.
1422 !%
1423 !% Note that the code will not check the orthonormality of the new states!
1424 !% Occupations are also fixed to the one of the previous states
1425 !%
1426 !% Each line provides the coefficients of the new states, in terms of
1427 !% the old ones. The coefficients are complex, but the imaginary part will be
1428 !% ignored for real wavefunctions.
1429 !% Note: This variable cannot be used when parallel in states.
1430 !%End
1431 if (parse_is_defined(namespace, trim(block_name))) then
1432 if (parse_block(namespace, trim(block_name), blk) == 0) then
1433 if (st%parallel_in_states) then
1434 call messages_not_implemented(trim(block_name) // " parallel in states", namespace=namespace)
1435 end if
1436 if (parse_block_n(blk) /= st%nst) then
1437 message(1) = "Number of rows in block " // trim(block_name) // " must equal number of states in this calculation."
1438 call messages_fatal(1, namespace=namespace)
1439 end if
1440 call states_elec_copy(stin, st, exclude_wfns = .true.)
1441 call states_elec_look_and_load(restart, namespace, space, stin, mesh, kpoints, fixed_occ = .true.)
1442
1443 ! FIXME: rotation matrix should be R_TYPE
1444 safe_allocate(rotation_matrix(1:stin%nst, 1:stin%nst))
1445 safe_allocate(psi(1:mesh%np, 1:st%d%dim))
1446
1447 rotation_matrix = diagonal_matrix(stin%nst, m_z1)
1448
1449 do ist = 1, st%nst
1450 ncols = parse_block_cols(blk, ist-1)
1451 if (ncols /= stin%nst) then
1452 write(message(1),'(a,i6,a,i6,3a,i6,a)') "Number of columns (", ncols, ") in row ", ist, " of block ", &
1453 trim(block_name), " must equal number of states (", stin%nst, ") read from gs restart."
1454 call messages_fatal(1, namespace=namespace)
1455 end if
1456 do jst = 1, stin%nst
1457 call parse_block_cmplx(blk, ist - 1, jst - 1, rotation_matrix(jst, ist))
1458 end do
1459 end do
1460
1461 call parse_block_end(blk)
1462
1463 do iqn = st%d%kpt%start, st%d%kpt%end
1464 if (states_are_real(st)) then
1465 call states_elec_rotate(stin, namespace, mesh, real(rotation_matrix, real64) , iqn)
1466 else
1467 call states_elec_rotate(stin, namespace, mesh, rotation_matrix, iqn)
1468 end if
1469
1470 do ist = st%st_start, st%st_end
1471 call states_elec_get_state(stin, mesh, ist, iqn, psi)
1472 call states_elec_set_state(st, mesh, ist, iqn, psi)
1473 end do
1474
1475 end do
1476
1477 safe_deallocate_a(rotation_matrix)
1478 safe_deallocate_a(psi)
1480 call states_elec_end(stin)
1481
1482 else
1483 call messages_variable_is_block(namespace, trim(block_name))
1484 end if
1485
1486 end if
1487
1488 pop_sub(states_elec_transform)
1489 end subroutine states_elec_transform
1490
1492
1493
1494!! Local Variables:
1495!! mode: f90
1496!! coding: utf-8
1497!! End:
ssize_t read(int __fd, void *__buf, size_t __nbytes) __attribute__((__access__(__write_only__
long int random(void)
Definition: getopt_f.c:1407
Copies a vector x, to a vector y.
Definition: lalg_basic.F90:186
This module handles the calculation mode.
integer, parameter, public p_strategy_max
integer, parameter, public p_strategy_domains
parallelization in domains
integer, parameter, public spinors
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_epsilon
Definition: global.F90:204
complex(real64), parameter, public m_z1
Definition: global.F90:199
subroutine, public zio_function_input(filename, namespace, space, mesh, ff, ierr, map)
Reads a mesh function from file filename, and puts it into ff. If the map argument is passed,...
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines various routines, operating on mesh functions.
subroutine, public zmf_normalize(mesh, dim, psi, norm)
Normalize a mesh function psi.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:921
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1114
character(len=512), private msg
Definition: messages.F90:166
subroutine, public messages_variable_is_block(namespace, name)
Definition: messages.F90:1070
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:538
subroutine, public print_date(str)
Definition: messages.F90:1006
subroutine, public messages_new_line()
Definition: messages.F90:1135
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
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
this module contains the low-level part of the output system
Definition: output_low.F90:115
character(len=max_path_len) function, public get_filename_with_spin(output, nspin, spin_index)
Returns the filame as output, or output-spX is spin polarized.
Definition: output_low.F90:232
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
subroutine, public restart_read(restart, iunit, lines, nlines, ierr)
Definition: restart.F90:935
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
subroutine, public smear_init(this, namespace, ispin, fixed_occ, integral_occs, kpoints)
Definition: smear.F90:185
pure logical function, public states_are_complex(st)
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
subroutine, public states_elec_end(st)
finalize the states_elec_t object
subroutine, public states_elec_allocate_wfns(st, mesh, wfs_type, skip, packed)
Allocates the KS wavefunctions defined within a states_elec_t structure.
subroutine, public states_elec_copy(stout, stin, exclude_wfns, exclude_eigenval, special)
make a (selective) copy of a states_elec_t object
subroutine, public states_elec_generate_random(st, mesh, kpoints, ist_start_, ist_end_, ikpt_start_, ikpt_end_, normalized)
randomize states
subroutine, public states_elec_look(restart, nik, dim, nst, ierr)
Reads the 'states' file in the restart directory, and finds out the nik, dim, and nst contained in it...
This module handles reading and writing restart information for the states_elec_t.
subroutine, public states_elec_read_user_def_orbitals(mesh, namespace, space, st)
the routine reads formulas for user-defined wavefunctions from the input file and fills the respectiv...
subroutine, public states_elec_load_frozen(restart, space, st, mesh, ierr)
subroutine, public states_elec_dump(restart, space, st, mesh, kpoints, ierr, iter, lr, st_start_writing, verbose)
subroutine, public states_elec_look_and_load(restart, namespace, space, st, mesh, kpoints, fixed_occ, is_complex, packed)
subroutine, public states_elec_transform(st, namespace, space, restart, mesh, kpoints, prefix)
subroutine, public states_elec_load(restart, namespace, space, st, mesh, kpoints, fixed_occ, ierr, iter, lr, lowest_missing, label, verbose, skip)
returns in ierr: <0 => Fatal error, or nothing read =0 => read all wavefunctions >0 => could only rea...
subroutine, public states_elec_load_spin(restart, st, ierr)
returns in ierr: <0 => Fatal error, or nothing read =0 => read all wavefunctions >0 => could only rea...
subroutine, public states_elec_dump_frozen(restart, space, st, mesh, ierr)
subroutine, public states_elec_load_rho(restart, space, st, mesh, ierr)
subroutine, public states_elec_dump_rho(restart, space, st, mesh, ierr, iter)
subroutine, public states_elec_dump_spin(restart, st, ierr)
subroutine, public conv_to_c_string(str)
converts to c string
Definition: string.F90:229
type(type_t), parameter, public type_cmplx
Definition: types.F90:134
type(type_t), parameter, public type_float
Definition: types.F90:133
subroutine fill_random()
logical function index_is_wrong()
Describes mesh distribution to nodes.
Definition: mesh.F90:186
The states_elec_t class contains all electronic wave functions.
int true(void)