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
41 use parser_oct_m
44 use smear_oct_m
45 use space_oct_m
50 use string_oct_m
51 use types_oct_m
52
53 implicit none
54
55
56 private
57 public :: &
69
70contains
71
72 ! ---------------------------------------------------------
73 subroutine states_elec_look_and_load(restart, namespace, space, st, mesh, kpoints, is_complex, packed)
74 type(restart_t), intent(in) :: restart
75 type(namespace_t), intent(in) :: namespace
76 class(space_t), intent(in) :: space
77 type(states_elec_t), target, intent(inout) :: st
78 class(mesh_t), intent(in) :: mesh
79 type(kpoints_t), intent(in) :: kpoints
80 logical, optional, intent(in) :: is_complex
81 logical, optional, intent(in) :: packed
82
83 integer :: nkpt, dim, nst, ierr
84 real(real64), allocatable :: new_occ(:,:)
85
87
88 !check how many wfs we have
89 call states_elec_look(restart, nkpt, dim, nst, ierr)
90 if (ierr /= 0) then
91 message(1) = "Unable to read states information."
92 call messages_fatal(1, namespace=namespace)
93 end if
94
95 if (st%parallel_in_states) then
96 message(1) = "Internal error: cannot use states_elec_look_and_load when parallel in states."
97 call messages_fatal(1, namespace=namespace)
98 end if
99
100 ! Resize st%occ, retaining information there
101 safe_allocate(new_occ(1:nst, 1:st%nik))
102 new_occ(:,:) = m_zero
103 new_occ(1:min(nst, st%nst),:) = st%occ(1:min(nst, st%nst),:)
104 safe_deallocate_a(st%occ)
105 call move_alloc(new_occ, st%occ)
106
107 ! FIXME: This wrong, one cannot just change the number of states
108 ! without updating the internal structures, in the case of parallel in states.
109 ! I think it is right now without state parallelism.
110 st%nst = nst
111 st%st_start = 1
112 st%st_end = nst
113 st%lnst = nst
114
115 safe_deallocate_a(st%node)
116 safe_allocate(st%node(1:st%nst))
117 st%node(:) = 0
118
119 safe_deallocate_a(st%eigenval)
120 safe_allocate(st%eigenval(1:st%nst, 1:st%nik))
121 st%eigenval = huge(st%eigenval)
122
123 if (present(is_complex)) then
124 if (is_complex) then
125 call states_elec_allocate_wfns(st, mesh, type_cmplx, packed=optional_default(packed, .false.))
126 else
127 call states_elec_allocate_wfns(st, mesh, type_float, packed=optional_default(packed, .false.))
128 end if
129 else
130 ! allow states_elec_allocate_wfns to decide for itself whether complex or real needed
131 call states_elec_allocate_wfns(st, mesh, packed=optional_default(packed, .false.))
132 end if
133
134 if (st%d%ispin == spinors) then
135 safe_allocate(st%spin(1:3, 1:st%nst, 1:st%nik))
136 st%spin = m_zero
137 end if
138
139 ! load wavefunctions
140 call states_elec_load(restart, namespace, space, st, mesh, kpoints, ierr)
141 if (ierr /= 0) then
142 message(1) = "Unable to read wavefunctions."
143 call messages_fatal(1, namespace=namespace)
144 end if
145
147 end subroutine states_elec_look_and_load
148
149
150 ! ---------------------------------------------------------
151 subroutine states_elec_dump(restart, space, st, mesh, kpoints, ierr, iter, lr, st_start_writing, verbose)
152 type(restart_t), intent(in) :: restart
153 class(space_t), intent(in) :: space
154 type(states_elec_t), intent(in) :: st
155 class(mesh_t), intent(in) :: mesh
156 type(kpoints_t), intent(in) :: kpoints
157 integer, intent(out) :: ierr
158 integer, optional, intent(in) :: iter
160 type(lr_t), optional, intent(in) :: lr
161 integer, optional, intent(in) :: st_start_writing
162 logical, optional, intent(in) :: verbose
163
164 integer :: iunit_wfns, iunit_occs, iunit_states
165 integer :: err, err2(2), ik, idir, ist, idim, itot
166 integer :: root(1:p_strategy_max)
167 character(len=MAX_PATH_LEN) :: filename
168 character(len=300) :: lines(3)
169 logical :: lr_wfns_are_associated, should_write, verbose_
170 real(real64) :: kpoint(space%dim)
171 real(real64), allocatable :: dpsi(:), rff_global(:)
172 complex(real64), allocatable :: zpsi(:), zff_global(:)
173
174 push_sub(states_elec_dump)
175
176 verbose_ = optional_default(verbose, .true.)
177
178 ierr = 0
179
180 if (restart_skip(restart)) then
181 pop_sub(states_elec_dump)
182 return
183 end if
184
185 if (verbose_) then
186 message(1) = "Info: Writing states."
187 call print_date(trim(message(1))//' ')
188 end if
189
190 call profiling_in("RESTART_WRITE")
191
192 if (present(lr)) then
193 lr_wfns_are_associated = (allocated(lr%ddl_psi) .and. states_are_real(st)) .or. &
194 (allocated(lr%zdl_psi) .and. states_are_complex(st))
195 assert(lr_wfns_are_associated)
196 end if
197
199
200
201 iunit_states = restart_open(restart, 'states')
202 write(lines(1), '(a20,1i10)') 'nst= ', st%nst
203 write(lines(2), '(a20,1i10)') 'dim= ', st%d%dim
204 write(lines(3), '(a20,1i10)') 'nik= ', st%nik
205 call restart_write(restart, iunit_states, lines, 3, err)
206 if (err /= 0) ierr = ierr + 1
207 call restart_close(restart, iunit_states)
208
209
210 iunit_wfns = restart_open(restart, 'wfns')
211 lines(1) = '# #k-point #st #dim filename'
212 if (states_are_real(st)) then
213 lines(2) = '%Real_Wavefunctions'
214 else
215 lines(2) = '%Complex_Wavefunctions'
216 end if
217 call restart_write(restart, iunit_wfns, lines, 2, err)
218 if (err /= 0) ierr = ierr + 2
219
220
221 iunit_occs = restart_open(restart, 'occs')
222 lines(1) = '# occupations | eigenvalue[a.u.] | Im(eigenvalue) [a.u.] | k-points | k-weights | filename | ik | ist | idim'
223 lines(2) = '%Occupations_Eigenvalues_K-Points'
224 call restart_write(restart, iunit_occs, lines, 2, err)
225 if (err /= 0) ierr = ierr + 4
226
227
228 if (states_are_real(st)) then
229 safe_allocate(dpsi(1:mesh%np))
230 safe_allocate(rff_global(1:mesh%np_global))
231 else
232 safe_allocate(zpsi(1:mesh%np))
233 safe_allocate(zff_global(1:mesh%np_global))
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 safe_deallocate_a(rff_global)
306 safe_deallocate_a(zff_global)
307
308 lines(1) = '%'
309 call restart_write(restart, iunit_occs, lines, 1, err)
310 if (err /= 0) ierr = ierr + 32
311 call restart_write(restart, iunit_wfns, lines, 1, err)
312 if (err /= 0) ierr = ierr + 64
313 if (present(iter)) then
314 write(lines(1),'(a,i7)') 'Iter = ', iter
315 call restart_write(restart, iunit_wfns, lines, 1, err)
316 if (err /= 0) ierr = ierr + 128
317 end if
318
319 call restart_close(restart, iunit_wfns)
320 call restart_close(restart, iunit_occs)
321
322 if (verbose_) then
323 message(1) = "Info: Finished writing states."
324 call print_date(trim(message(1))//' ')
325 end if
326
328
329 call profiling_out("RESTART_WRITE")
330 pop_sub(states_elec_dump)
331 return
332 end subroutine states_elec_dump
333
334
335 ! ---------------------------------------------------------
340 subroutine states_elec_load(restart, namespace, space, st, mesh, kpoints, ierr, iter, lr, lowest_missing, label, verbose, skip)
341 type(restart_t), intent(in) :: restart
342 type(namespace_t), intent(in) :: namespace
343 class(space_t), intent(in) :: space
344 type(states_elec_t), intent(inout) :: st
345 class(mesh_t), intent(in) :: mesh
346 type(kpoints_t), intent(in) :: kpoints
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(:), zpsil(:)
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 (st%restart_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
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(zpsil)
672 safe_deallocate_a(restart_file)
673 safe_deallocate_a(restart_file_present)
674
675 if (mpi_grp_is_root(mpi_world) .and. verbose_) then
676 call messages_new_line()
677 end if
678
679 if (st%parallel_in_states .or. st%d%kpt%parallel) then
680 iread_tmp = iread
681 call st%st_kpt_mpi_grp%allreduce(iread_tmp, iread, 1, mpi_integer, mpi_sum)
682 end if
683
684 if (st%d%kpt%parallel) then
685 ! make sure all tasks know lowest_missing over all k-points
686 if (present(lowest_missing)) then
687 safe_allocate(lowest_missing_tmp(1:st%d%dim, 1:st%nik))
688 lowest_missing_tmp = lowest_missing
689 call st%d%kpt%mpi_grp%allreduce(lowest_missing_tmp(1,1), lowest_missing(1,1), st%d%dim*st%nik, &
690 mpi_integer, mpi_min)
691 safe_deallocate_a(lowest_missing_tmp)
692 end if
693 end if
694
695 if (st%restart_fixed_occ .and. iread == st%nst * st%nik * st%d%dim) then
696 ! reset to overwrite whatever smearing may have been set earlier
697 ! only do this if all states could be loaded from restart files
698 call smear_init(st%smear, namespace, st%d%ispin, fixed_occ = .true., integral_occs = integral_occs, kpoints = kpoints)
699 end if
700
701 if (.not. present(lr) .and. .not. present(skip)) call fill_random()
702 ! it is better to initialize lr wfns to zero
703
704 safe_deallocate_a(filled)
705
706 if (ierr == 0 .and. iread /= st%nst * st%nik * st%d%dim) then
707 if (iread > 0) then
708 ierr = iread
709 else
710 ierr = -1
711 end if
712 ! otherwise ierr = 0 would mean either all was read correctly, or nothing at all was read!
713
714 if (.not. present(lr)) then
715 write(str, '(a,i5)') 'Reading states.'
716 else
717 write(str, '(a,i5)') 'Reading states information for linear response.'
718 end if
719 call messages_print_with_emphasis(msg=trim(str), namespace=namespace)
720 if (.not. present(skip)) then
721 write(message(1),'(a,i6,a,i6,a)') 'Only ', iread,' files out of ', &
722 st%nst * st%nik * st%d%dim, ' could be read.'
723 else
724 write(message(1),'(a,i6,a,i6,a)') 'Only ', iread,' files out of ', &
725 st%nst * st%nik * st%d%dim, ' were loaded.'
726 ierr = 0
727 end if
728 call messages_info(1, namespace=namespace)
729 call messages_print_with_emphasis(namespace=namespace)
730 end if
731
732 message(1) = 'Info: States reading done.'
733 if (verbose_) call print_date(trim(message(1))//' ')
734
735 call profiling_out("RESTART_READ")
736 pop_sub(states_elec_load)
737
738 contains
739
740 ! ---------------------------------------------------------
741 subroutine fill_random() !< put random function in orbitals that could not be read.
743
744 do ik = st%d%kpt%start, st%d%kpt%end
745
746 do ist = st%st_start, st%st_end
747 do idim = 1, st%d%dim
748 if (filled(idim, ist, ik)) cycle
749
750 call states_elec_generate_random(st, mesh, kpoints, ist, ist, ik, ik)
751 end do
752 end do
753 end do
754
756 end subroutine fill_random
757
758 ! ---------------------------------------------------------
759
760 logical function index_is_wrong() !< .true. if the index (idim, ist, ik) is not present in st structure...
762
763 if (idim > st%d%dim .or. idim < 1 .or. &
764 ist > st%nst .or. ist < 1 .or. &
765 ik > st%nik .or. ik < 1) then
767 else
768 index_is_wrong = .false.
769 end if
770
772 end function index_is_wrong
773
774 end subroutine states_elec_load
775
776
777 subroutine states_elec_dump_rho(restart, space, st, mesh, ierr, iter)
778 type(restart_t), intent(in) :: restart
779 class(space_t), intent(in) :: space
780 type(states_elec_t), intent(in) :: st
781 class(mesh_t), intent(in) :: mesh
782 integer, intent(out) :: ierr
783 integer, optional, intent(in) :: iter
784
785 integer :: iunit, isp, err, err2(2)
786 character(len=80) :: filename
787 character(len=300) :: lines(2)
788
789 push_sub(states_elec_dump_rho)
790
791 ierr = 0
792
793 if (restart_skip(restart)) then
794 pop_sub(states_elec_dump_rho)
795 return
796 end if
797
798 if (debug%info) then
799 message(1) = "Debug: Writing density restart."
800 call messages_info(1)
801 end if
802
804
805 !write the densities
806 iunit = restart_open(restart, 'density')
807 lines(1) = '# #spin #nspin filename'
808 lines(2) = '%densities'
809 call restart_write(restart, iunit, lines, 2, err)
810 if (err /= 0) ierr = ierr + 1
811
812 err2 = 0
813 do isp = 1, st%d%nspin
814 if (st%d%nspin == 1) then
815 write(filename, fmt='(a)') 'density'
816 else
817 write(filename, fmt='(a,i1)') 'density-sp', isp
818 end if
819 write(lines(1), '(i8,a,i8,a)') isp, ' | ', st%d%nspin, ' | "'//trim(adjustl(filename))//'"'
820 call restart_write(restart, iunit, lines, 1, err)
821 if (err /= 0) err2(1) = err2(1) + 1
822
823 call restart_write_mesh_function(restart, space, filename, mesh, st%rho(:,isp), err)
824 if (err /= 0) err2(2) = err2(2) + 1
825
826 end do
827 if (err2(1) /= 0) ierr = ierr + 2
828 if (err2(2) /= 0) ierr = ierr + 4
829
830 lines(1) = '%'
831 call restart_write(restart, iunit, lines, 1, err)
832 if (err /= 0) ierr = ierr + 8
833 if (present(iter)) then
834 write(lines(1),'(a,i7)') 'Iter = ', iter
835 call restart_write(restart, iunit, lines, 1, err)
836 if (err /= 0) ierr = ierr + 16
837 end if
838 call restart_close(restart, iunit)
839
841
842 if (debug%info) then
843 message(1) = "Debug: Writing density restart done."
844 call messages_info(1)
845 end if
846
847 pop_sub(states_elec_dump_rho)
848 end subroutine states_elec_dump_rho
849
850
851 ! ---------------------------------------------------------
852 subroutine states_elec_load_rho(restart, space, st, mesh, ierr)
853 type(restart_t), intent(in) :: restart
854 class(space_t), intent(in) :: space
855 type(states_elec_t), intent(inout) :: st
856 class(mesh_t), intent(in) :: mesh
857 integer, intent(out) :: ierr
858
859 integer :: err, err2, isp
860 character(len=12) :: filename
861
862 push_sub(states_elec_load_rho)
863
864 ierr = 0
865
866 if (restart_skip(restart)) then
867 ierr = -1
868 pop_sub(states_elec_load_rho)
869 return
870 end if
871
872 if (debug%info) then
873 message(1) = "Debug: Reading density restart."
874 call messages_info(1)
875 end if
876
877 ! skip for now, since we know what the files are going to be called
878 !read the densities
879! iunit_rho = io_open(trim(dir)//'/density', restart%namespace, action='write')
880! call iopar_read(st%dom_st_kpt_mpi_grp, iunit_rho, line, err)
881! call iopar_read(st%dom_st_kpt_mpi_grp, iunit_rho, line, err)
882! we could read the iteration 'iter' too, not sure if that is useful.
883
884 err2 = 0
885 do isp = 1, st%d%nspin
886 if (st%d%nspin == 1) then
887 write(filename, fmt='(a)') 'density'
888 else
889 write(filename, fmt='(a,i1)') 'density-sp', isp
890 end if
891! if (mpi_grp_is_root(st%dom_st_kpt_mpi_grp)) then
892! read(iunit_rho, '(i8,a,i8,a)') isp, ' | ', st%d%nspin, ' | "'//trim(adjustl(filename))//'"'
893! end if
894 call restart_read_mesh_function(restart, space, filename, mesh, st%rho(:,isp), err)
895 if (err /= 0) err2 = err2 + 1
896
897 end do
898 if (err2 /= 0) ierr = ierr + 1
899
900 if (debug%info) then
901 message(1) = "Debug: Reading density restart done."
902 call messages_info(1)
903 end if
904
905 pop_sub(states_elec_load_rho)
906 end subroutine states_elec_load_rho
907
908 subroutine states_elec_dump_frozen(restart, space, st, mesh, ierr)
909 type(restart_t), intent(in) :: restart
910 class(space_t), intent(in) :: space
911 type(states_elec_t), intent(in) :: st
912 class(mesh_t), intent(in) :: mesh
913 integer, intent(out) :: ierr
914
915 integer :: isp, err, err2(2), idir
916 character(len=80) :: filename
917
919
920 ierr = 0
921
922 assert(allocated(st%frozen_rho))
923
924 if (restart_skip(restart)) then
926 return
927 end if
928
929 if (debug%info) then
930 message(1) = "Debug: Writing frozen densities restart."
931 call messages_info(1)
932 end if
933
935
936 err2 = 0
937 do isp = 1, st%d%nspin
938 if (st%d%nspin == 1) then
939 write(filename, fmt='(a)') 'frozen_rho'
940 else
941 write(filename, fmt='(a,i1)') 'frozen_rho-sp', isp
942 end if
943
944 call restart_write_mesh_function(restart, space, filename, mesh, st%frozen_rho(:,isp), err)
945 if (err /= 0) err2(2) = err2(2) + 1
946
947 if (allocated(st%frozen_tau)) then
948 if (st%d%nspin == 1) then
949 write(filename, fmt='(a)') 'frozen_tau'
950 else
951 write(filename, fmt='(a,i1)') 'frozen_tau-sp', isp
952 end if
953 call restart_write_mesh_function(restart, space, filename, mesh, st%frozen_tau(:,isp), err)
954 if (err /= 0) err2 = err2 + 1
955 end if
956
957 if (allocated(st%frozen_gdens)) then
958 do idir = 1, space%dim
959 if (st%d%nspin == 1) then
960 write(filename, fmt='(a,i1)') 'frozen_gdens-dir', idir
961 else
962 write(filename, fmt='(a,i1,a,i1)') 'frozen_tau-dir', idir, '-', isp
963 end if
964 call restart_write_mesh_function(restart, space, filename, mesh, st%frozen_gdens(:,idir,isp), err)
965 if (err /= 0) err2 = err2 + 1
966 end do
967 end if
968
969 if (allocated(st%frozen_ldens)) then
970 if (st%d%nspin == 1) then
971 write(filename, fmt='(a)') 'frozen_ldens'
972 else
973 write(filename, fmt='(a,i1)') 'frozen_ldens-sp', isp
974 end if
975 call restart_write_mesh_function(restart, space, filename, mesh, st%frozen_ldens(:,isp), err)
976 if (err /= 0) err2 = err2 + 1
977 end if
978
979 end do
980 if (err2(1) /= 0) ierr = ierr + 2
981 if (err2(2) /= 0) ierr = ierr + 4
982
984
985 if (debug%info) then
986 message(1) = "Debug: Writing frozen densities restart done."
987 call messages_info(1)
988 end if
989
991 end subroutine states_elec_dump_frozen
992
993
994 ! ---------------------------------------------------------
995 subroutine states_elec_load_frozen(restart, space, st, mesh, ierr)
996 type(restart_t), intent(in) :: restart
997 class(space_t), intent(in) :: space
998 type(states_elec_t), intent(inout) :: st
999 class(mesh_t), intent(in) :: mesh
1000 integer, intent(out) :: ierr
1002 integer :: err, err2, isp, idir
1003 character(len=80) :: filename
1004
1005 push_sub(states_elec_load_frozen)
1006
1007 assert(allocated(st%frozen_rho))
1008
1009 ierr = 0
1010
1011 if (restart_skip(restart)) then
1012 ierr = -1
1014 return
1015 end if
1016
1017 if (debug%info) then
1018 message(1) = "Debug: Reading densities restart."
1019 call messages_info(1)
1020 end if
1021
1022 err2 = 0
1023 do isp = 1, st%d%nspin
1024 if (st%d%nspin == 1) then
1025 write(filename, fmt='(a)') 'frozen_rho'
1026 else
1027 write(filename, fmt='(a,i1)') 'frozen_rho-sp', isp
1028 end if
1029 call restart_read_mesh_function(restart, space, filename, mesh, st%frozen_rho(:,isp), err)
1030 if (err /= 0) err2 = err2 + 1
1031
1032 if (allocated(st%frozen_tau)) then
1033 if (st%d%nspin == 1) then
1034 write(filename, fmt='(a)') 'frozen_tau'
1035 else
1036 write(filename, fmt='(a,i1)') 'frozen_tau-sp', isp
1037 end if
1038 call restart_read_mesh_function(restart, space, filename, mesh, st%frozen_tau(:,isp), err)
1039 if (err /= 0) err2 = err2 + 1
1040 end if
1041
1042 if (allocated(st%frozen_gdens)) then
1043 do idir = 1, space%dim
1044 if (st%d%nspin == 1) then
1045 write(filename, fmt='(a,i1)') 'frozen_gdens-dir', idir
1046 else
1047 write(filename, fmt='(a,i1,a,i1)') 'frozen_tau-dir', idir, '-', isp
1048 end if
1049 call restart_read_mesh_function(restart, space, filename, mesh, st%frozen_gdens(:,idir,isp), err)
1050 if (err /= 0) err2 = err2 + 1
1051 end do
1052 end if
1053
1054 if (allocated(st%frozen_ldens)) then
1055 if (st%d%nspin == 1) then
1056 write(filename, fmt='(a)') 'frozen_ldens'
1057 else
1058 write(filename, fmt='(a,i1)') 'frozen_ldens-sp', isp
1059 end if
1060 call restart_read_mesh_function(restart, space, filename, mesh, st%frozen_ldens(:,isp), err)
1061 if (err /= 0) err2 = err2 + 1
1062 end if
1063
1064 end do
1065 if (err2 /= 0) ierr = ierr + 1
1066
1067 if (debug%info) then
1068 message(1) = "Debug: Reading frozen densities restart done."
1069 call messages_info(1)
1070 end if
1071
1073 end subroutine states_elec_load_frozen
1074
1075
1076 ! ---------------------------------------------------------
1079 subroutine states_elec_read_user_def_orbitals(mesh, namespace, space, st)
1080 class(mesh_t), intent(in) :: mesh
1081 type(namespace_t), intent(in) :: namespace
1082 class(space_t), intent(in) :: space
1083 type(states_elec_t), intent(inout) :: st
1084
1085 type(block_t) :: blk
1086 integer :: ip, id, is, ik, nstates, state_from, ierr, ncols
1087 integer :: ib, idim, inst, inik, normalize
1088 real(real64) :: xx(space%dim), rr, psi_re, psi_im
1089 character(len=150) :: filename
1090 complex(real64), allocatable :: zpsi(:, :)
1091
1092 integer, parameter :: &
1093 state_from_formula = 1, &
1094 state_from_file = -10010, &
1095 normalize_yes = 1, &
1096 normalize_no = 0
1097
1099
1100 !%Variable UserDefinedStates
1101 !%Type block
1102 !%Section States
1103 !%Description
1104 !% Instead of using the ground state as initial state for
1105 !% time-propagations it might be interesting in some cases
1106 !% to specify alternate states. Like with user-defined
1107 !% potentials, this block allows you to specify formulas for
1108 !% the orbitals at <i>t</i>=0.
1109 !%
1110 !% Example:
1111 !%
1112 !% <tt>%UserDefinedStates
1113 !% <br>&nbsp;&nbsp; 1 | 1 | 1 | formula | "exp(-r^2)*exp(-i*0.2*x)" | normalize_yes
1114 !% <br>%</tt>
1115 !%
1116 !% The first column specifies the component of the spinor,
1117 !% the second column the number of the state and the third
1118 !% contains <i>k</i>-point and spin quantum numbers. Column four
1119 !% indicates that column five should be interpreted as a formula
1120 !% for the corresponding orbital.
1121 !%
1122 !% Alternatively, if column four states <tt>file</tt> the state will
1123 !% be read from the file given in column five.
1124 !%
1125 !% <tt>%UserDefinedStates
1126 !% <br>&nbsp;&nbsp; 1 | 1 | 1 | file | "/path/to/file" | normalize_no
1127 !% <br>%</tt>
1128 !%
1129 !% Octopus reads first the ground-state orbitals from
1130 !% the <tt>restart/gs</tt> directory. Only the states that are
1131 !% specified in the above block will be overwritten with
1132 !% the given analytic expression for the orbital.
1133 !%
1134 !% The sixth (optional) column indicates whether <tt>Octopus</tt> should renormalize
1135 !% the orbital. The default (no sixth column given) is to renormalize.
1136 !%
1137 !%Option file -10010
1138 !% Read initial orbital from file.
1139 !% Accepted file formats, detected by extension: obf, ncdf and csv (real only).
1140 !%Option formula 1
1141 !% Calculate initial orbital by given analytic expression.
1142 !%Option normalize_yes 1
1143 !% Normalize orbitals (default).
1144 !%Option normalize_no 0
1145 !% Do not normalize orbitals.
1146 !%End
1147 if (parse_block(namespace, 'UserDefinedStates', blk) == 0) then
1148
1149 call messages_print_with_emphasis(msg=trim('Substitution of orbitals'), namespace=namespace)
1150
1151 ! find out how many lines (i.e. states) the block has
1152 nstates = parse_block_n(blk)
1153
1154 safe_allocate(zpsi(1:mesh%np, 1:st%d%dim))
1155
1156 ! read all lines
1157 do ib = 1, nstates
1158 ! Check that number of columns is five or six.
1159 ncols = parse_block_cols(blk, ib - 1)
1160 if (ncols < 5 .or. ncols > 6) then
1161 message(1) = 'Each line in the UserDefinedStates block must have'
1162 message(2) = 'five or six columns.'
1163 call messages_fatal(2, namespace=namespace)
1164 end if
1165
1166 call parse_block_integer(blk, ib - 1, 0, idim)
1167 call parse_block_integer(blk, ib - 1, 1, inst)
1168 call parse_block_integer(blk, ib - 1, 2, inik)
1169
1170 ! Calculate from expression or read from file?
1171 call parse_block_integer(blk, ib - 1, 3, state_from)
1173 ! loop over all states
1174 do id = 1, st%d%dim
1175 do is = 1, st%nst
1176 do ik = 1, st%nik
1177
1178 if (.not.(id == idim .and. is == inst .and. ik == inik )) cycle
1179
1180 ! We first parse the value, such that the parsing and stdout looks correct
1181 ! However, if rank 0 is not responsible for the parsing, the variables in the formula
1182 ! will still appear as not parsed (pure k-point/state parallelization case).
1183 ! This could be solved by a dummy call to parse_expression
1184 select case (state_from)
1185 case (state_from_formula)
1186 ! parse formula string
1187 call parse_block_string( &
1188 blk, ib - 1, 4, st%user_def_states(id, is, ik))
1189
1190 write(message(1), '(a,3i5)') 'Substituting state of orbital with k, ist, dim = ', ik, is, id
1191 write(message(2), '(2a)') ' with the expression:'
1192 write(message(3), '(2a)') ' ',trim(st%user_def_states(id, is, ik))
1193 call messages_info(3, namespace=namespace)
1194
1195 ! convert to C string
1196 call conv_to_c_string(st%user_def_states(id, is, ik))
1197
1198 case (state_from_file)
1199 ! The input format can be coded in column four now. As it is
1200 ! not used now, we just say "file".
1201 ! Read the filename.
1202 call parse_block_string(blk, ib - 1, 4, filename)
1203
1204 write(message(1), '(a,3i5)') 'Substituting state of orbital with k, ist, dim = ', ik, is, id
1205 write(message(2), '(2a)') ' with data from file:'
1206 write(message(3), '(2a)') ' ',trim(filename)
1207 call messages_info(3, namespace=namespace)
1208 end select
1209
1210
1211 ! does the block entry match and is this node responsible?
1212 if (.not.(st%st_start <= is .and. st%st_end >= is &
1213 .and. st%d%kpt%start <= ik .and. st%d%kpt%end >= ik)) cycle
1214
1215 select case (state_from)
1216
1217 case (state_from_formula)
1218 ! fill states with user-defined formulas
1219 do ip = 1, mesh%np
1220 xx = mesh%x(ip, :)
1221 rr = norm2(xx)
1222
1223 ! parse user-defined expressions
1224 call parse_expression(psi_re, psi_im, space%dim, xx, rr, m_zero, st%user_def_states(id, is, ik))
1225 ! fill state
1226 zpsi(ip, 1) = psi_re + m_zi * psi_im
1227 end do
1228
1229 case (state_from_file)
1230 ! finally read the state
1231 call zio_function_input(filename, namespace, space, mesh, zpsi(:, 1), ierr)
1232 if (ierr > 0) then
1233 message(1) = 'Could not read the file!'
1234 write(message(2),'(a,i1)') 'Error code: ', ierr
1235 call messages_fatal(2, namespace=namespace)
1236 end if
1237
1238 case default
1239 message(1) = 'Wrong entry in UserDefinedStates, column 4.'
1240 message(2) = 'You may state "formula" or "file" here.'
1241 call messages_fatal(2, namespace=namespace)
1242 end select
1243
1244 ! normalize orbital
1245 if (parse_block_cols(blk, ib - 1) == 6) then
1246 call parse_block_integer(blk, ib - 1, 5, normalize)
1247 else
1248 normalize = 1
1249 end if
1250 select case (normalize)
1251 case (normalize_no)
1252 call states_elec_set_state(st, mesh, id, is, ik, zpsi(:, 1))
1253 case (normalize_yes)
1254 assert(st%d%dim == 1)
1255 call zmf_normalize(mesh, st%d%dim, zpsi)
1256 call states_elec_set_state(st, mesh, is, ik, zpsi)
1257 case default
1258 message(1) = 'The sixth column in UserDefinedStates may either be'
1259 message(2) = '"normalize_yes" or "normalize_no"'
1260 call messages_fatal(2, namespace=namespace)
1261 end select
1262
1263 end do
1264 end do
1265 end do
1266
1267 end do
1268
1269 safe_deallocate_a(zpsi)
1270
1271 call parse_block_end(blk)
1272 call messages_print_with_emphasis(namespace=namespace)
1273
1274 else
1275 call messages_variable_is_block(namespace, 'UserDefinedStates')
1276 end if
1277
1280
1281 ! ---------------------------------------------------------
1282 ! This is needed as for the generalized Bloch theorem we need to label
1283 ! the states from the expectation value of Sz computed from the GS.
1284 subroutine states_elec_dump_spin(restart, st, ierr)
1285 type(restart_t), intent(in) :: restart
1286 type(states_elec_t), intent(in) :: st
1287 integer, intent(out) :: ierr
1288
1289 integer :: iunit_spin
1290 integer :: err, err2(2), ik, ist
1291 character(len=300) :: lines(3)
1292
1293 push_sub(states_elec_dump_spin)
1294
1295 ierr = 0
1296
1297 if (restart_skip(restart)) then
1298 pop_sub(states_elec_dump_spin)
1299 return
1300 end if
1301
1302 call profiling_in("RESTART_WRITE_SPIN")
1303
1305
1306 iunit_spin = restart_open(restart, 'spin')
1307 lines(1) = '# #k-point #st #spin(x) spin(y) spin(z)'
1308 call restart_write(restart, iunit_spin, lines, 1, err)
1309 if (err /= 0) ierr = ierr + 1
1310
1311 err2 = 0
1312 do ik = 1, st%nik
1313 do ist = 1, st%nst
1314 write(lines(1), '(i8,a,i8,3(a,f18.12))') ik, ' | ', ist, ' | ', &
1315 st%spin(1,ist,ik), ' | ', st%spin(2,ist,ik),' | ', st%spin(3,ist,ik)
1316 call restart_write(restart, iunit_spin, lines, 1, err)
1317 if (err /= 0) err2(1) = err2(1) + 1
1318 end do ! st%nst
1319 end do ! st%nik
1320 lines(1) = '%'
1321 call restart_write(restart, iunit_spin, lines, 1, err)
1322 if (err2(1) /= 0) ierr = ierr + 8
1323 if (err2(2) /= 0) ierr = ierr + 16
1324
1325 call restart_close(restart, iunit_spin)
1326
1328
1329 call profiling_out("RESTART_WRITE_SPIN")
1330 pop_sub(states_elec_dump_spin)
1331 end subroutine states_elec_dump_spin
1332
1333
1334
1335
1336 ! ---------------------------------------------------------
1341 subroutine states_elec_load_spin(restart, st, ierr)
1342 type(restart_t), intent(in) :: restart
1343 type(states_elec_t), intent(inout) :: st
1344 integer, intent(out) :: ierr
1345
1346 integer :: spin_file, err, ik, ist
1347 character(len=256) :: lines(3)
1348 real(real64) :: spin(3)
1349 character(len=1) :: char
1350
1351
1352 push_sub(states_elec_load_spin)
1353
1354 ierr = 0
1355
1356 if (restart_skip(restart)) then
1357 ierr = -1
1358 pop_sub(states_elec_load_spin)
1359 return
1360 end if
1361
1362 call profiling_in("RESTART_READ_SPIN")
1363
1364 ! open files to read
1365 spin_file = restart_open(restart, 'spin')
1366 ! Skip two lines.
1367 call restart_read(restart, spin_file, lines, 1, err)
1368 if (err /= 0) ierr = ierr - 2**7
1369
1370 ! If any error occured up to this point then it is not worth continuing,
1371 ! as there something fundamentally wrong with the restart files
1372 if (ierr /= 0) then
1373 call restart_close(restart, spin_file)
1374 call profiling_out("RESTART_READ_SPIN")
1375 pop_sub(states_elec_load_spin)
1376 return
1377 end if
1378
1379 ! Next we read the list of states from the files.
1380 ! Errors in reading the information of a specific state from the files are ignored
1381 ! at this point, because later we will skip reading the wavefunction of that state.
1382 do
1383 call restart_read(restart, spin_file, lines, 1, err)
1384 read(lines(1), '(a)') char
1385 if (char == '%') then
1386 !We reached the end of the file
1387 exit
1388 else
1389 read(lines(1), *) ik, char, ist, char, spin(1), char, spin(2), char, spin(3)
1390 end if
1391
1392 if (err /= 0) cycle
1393
1394 st%spin(1:3, ist, ik) = spin(1:3)
1395 end do
1396
1397 call restart_close(restart, spin_file)
1398
1399 call profiling_out("RESTART_READ_SPIN")
1400 pop_sub(states_elec_load_spin)
1401 end subroutine states_elec_load_spin
1402
1403 ! ---------------------------------------------------------
1404 subroutine states_elec_transform(st, namespace, space, restart, mesh, kpoints, prefix)
1405 type(states_elec_t), intent(inout) :: st
1406 type(namespace_t), intent(in) :: namespace
1407 class(space_t), intent(in) :: space
1408 type(restart_t), intent(inout) :: restart
1409 class(mesh_t), intent(in) :: mesh
1410 type(kpoints_t), intent(in) :: kpoints
1411 character(len=*), optional, intent(in) :: prefix
1412
1413 type(states_elec_t) :: stin
1414 type(block_t) :: blk
1415 complex(real64), allocatable :: rotation_matrix(:,:), psi(:, :)
1416 integer :: ist, jst, ncols, iqn
1417 character(len=256) :: block_name
1418
1419 push_sub(states_elec_transform)
1420
1421 block_name = trim(optional_default(prefix, "")) // "TransformStates"
1422
1423 !%Variable TransformStates
1424 !%Type block
1425 !%Default no
1426 !%Section States
1427 !%Description
1428 !% Before starting the <tt>td</tt> calculation, the initial states (that are
1429 !% read from the <tt>restart/gs</tt> directory, which should have been
1430 !% generated in a previous ground-state calculation) can be "transformed"
1431 !% among themselves. The block <tt>TransformStates</tt> gives the transformation matrix
1432 !% to be used. The number of rows and columns of the matrix should equal the number
1433 !% of the states present in the time-dependent calculation (the independent
1434 !% spin and <i>k</i>-point subspaces are all transformed equally); the number of
1435 !% columns should be equal to the number of states present in the
1436 !% <tt>restart/gs</tt> directory. This number may be different: for example,
1437 !% one could have run previously in <tt>unocc</tt> mode in order to obtain unoccupied
1438 !% Kohn-Sham states, and therefore <tt>restart/gs</tt> will contain more states.
1439 !% These states can be used in the transformation.
1440 !%
1441 !% Note that the code will not check the orthonormality of the new states!
1442 !%
1443 !% Each line provides the coefficients of the new states, in terms of
1444 !% the old ones. The coefficients are complex, but the imaginary part will be
1445 !% ignored for real wavefunctions.
1446 !% Note: This variable cannot be used when parallel in states.
1447 !%End
1448 if (parse_is_defined(namespace, trim(block_name))) then
1449 if (parse_block(namespace, trim(block_name), blk) == 0) then
1450 if (st%parallel_in_states) then
1451 call messages_not_implemented(trim(block_name) // " parallel in states", namespace=namespace)
1452 end if
1453 if (parse_block_n(blk) /= st%nst) then
1454 message(1) = "Number of rows in block " // trim(block_name) // " must equal number of states in this calculation."
1455 call messages_fatal(1, namespace=namespace)
1456 end if
1457 call states_elec_copy(stin, st, exclude_wfns = .true.)
1458 call states_elec_look_and_load(restart, namespace, space, stin, mesh, kpoints)
1459
1460 ! FIXME: rotation matrix should be R_TYPE
1461 safe_allocate(rotation_matrix(1:stin%nst, 1:stin%nst))
1462 safe_allocate(psi(1:mesh%np, 1:st%d%dim))
1463
1464 rotation_matrix = diagonal_matrix(stin%nst, m_z1)
1465
1466 do ist = 1, st%nst
1467 ncols = parse_block_cols(blk, ist-1)
1468 if (ncols /= stin%nst) then
1469 write(message(1),'(a,i6,a,i6,3a,i6,a)') "Number of columns (", ncols, ") in row ", ist, " of block ", &
1470 trim(block_name), " must equal number of states (", stin%nst, ") read from gs restart."
1471 call messages_fatal(1, namespace=namespace)
1472 end if
1473 do jst = 1, stin%nst
1474 call parse_block_cmplx(blk, ist - 1, jst - 1, rotation_matrix(jst, ist))
1475 end do
1476 end do
1477
1478 call parse_block_end(blk)
1479
1480 do iqn = st%d%kpt%start, st%d%kpt%end
1481 if (states_are_real(st)) then
1482 call states_elec_rotate(stin, namespace, mesh, real(rotation_matrix, real64) , iqn)
1483 else
1484 call states_elec_rotate(stin, namespace, mesh, rotation_matrix, iqn)
1485 end if
1486
1487 do ist = st%st_start, st%st_end
1488 call states_elec_get_state(stin, mesh, ist, iqn, psi)
1489 call states_elec_set_state(st, mesh, ist, iqn, psi)
1490 end do
1491
1492 end do
1493
1494 safe_deallocate_a(rotation_matrix)
1495 safe_deallocate_a(psi)
1496
1498
1499 else
1500 call messages_variable_is_block(namespace, trim(block_name))
1501 end if
1502
1503 end if
1504
1505 pop_sub(states_elec_transform)
1506 end subroutine states_elec_transform
1507
1509
1510
1511!! Local Variables:
1512!! mode: f90
1513!! coding: utf-8
1514!! End:
long int random(void)
Definition: getopt_f.c:1407
Copies a vector x, to a vector y.
Definition: lalg_basic.F90:185
This module handles the calculation mode.
integer, parameter, public p_strategy_max
integer, parameter, public p_strategy_domains
parallelization in domains
type(debug_t), save, public debug
Definition: debug.F90:154
integer, parameter, public spinors
real(real64), parameter, public m_zero
Definition: global.F90:187
complex(real64), parameter, public m_zi
Definition: global.F90:201
real(real64), parameter, public m_epsilon
Definition: global.F90:203
complex(real64), parameter, public m_z1
Definition: global.F90:198
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:930
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1125
character(len=512), private msg
Definition: messages.F90:165
subroutine, public messages_variable_is_block(namespace, name)
Definition: messages.F90:1081
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
subroutine, public print_date(str)
Definition: messages.F90:1017
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
subroutine, public messages_new_line()
Definition: messages.F90:1146
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
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
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:933
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
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, 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, 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:252
type(type_t), public type_float
Definition: types.F90:133
type(type_t), public type_cmplx
Definition: types.F90:134
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)