56 type(electrons_t),
pointer :: sys
57 type(states_elec_t) :: st
58 complex(real64),
allocatable :: hmss(:,:), psi(:,:,:), hpsi(:,:,:), temp_state1(:,:)
59 complex(real64),
allocatable :: HFloquet(:,:,:), HFloq_eff(:,:), temp(:,:)
60 real(real64),
allocatable :: eigenval(:), bands(:,:)
61 character(len=80) :: filename
62 integer :: it, nT, ik, ist, in, im, file, idim, nik, ik_count
63 integer :: Forder, Fdim, m0, n0, n1, nst, ii, jj, lim_nst
64 real(real64) :: dt, Tcycle,omega
65 logical :: downfolding = .false.
66 class(mesh_t),
pointer :: mesh
67 type(restart_t) :: restart
94 sys%ext_partners, st, time=
m_zero)
104 message(1) =
'Unable to read ground-state wavefunctions.'
122 safe_deallocate_p(sys)
139 assert(sys%gr%np == sys%gr%np_global)
145 message(1) =
"Please give a non-zero value for TDFloquetFrequency"
154 dt = tcycle/real(nt, real64)
157 if (forder.ge.0)
then
160 fdim = 2 * forder + 1
162 message(1) =
'Floquet-Hamiltonian is downfolded'
169 dt = tcycle/real(nt, real64)
184 safe_allocate(hmss(1:nst,1:nst))
185 safe_allocate( psi(1:nst,1:st%d%dim,1:mesh%np))
186 safe_allocate(hpsi(1:nst,1:st%d%dim,1:mesh%np))
187 safe_allocate(temp_state1(1:mesh%np,1:st%d%dim))
193 nik = sys%kpoints%nik_skip
201 safe_allocate(hfloquet(1:nik,1:nst*fdim, 1:nst*fdim))
202 hfloquet(1:nik,1:nst*fdim, 1:nst*fdim) =
m_zero
215 do ik = sys%kpoints%reduced%npoints-nik+1, sys%kpoints%reduced%npoints
216 ik_count = ik_count + 1
218 psi(1:nst, 1:st%d%dim, 1:mesh%np)=
m_zero
219 hpsi(1:nst, 1:st%d%dim, 1:mesh%np)=
m_zero
221 do ist = st%st_start, st%st_end
224 do idim = 1, st%d%dim
225 psi(ist, idim, 1:mesh%np) = temp_state1(1:mesh%np, idim)
228 do idim = 1, st%d%dim
229 hpsi(ist, idim, 1:mesh%np) = temp_state1(1:mesh%np, idim)
235 assert(mesh%np_global*st%d%dim < huge(0_int32))
236 hmss(1:nst, 1:nst) =
m_zero
241 i8_to_i4(mesh%np_global*st%d%dim), &
242 cmplx(mesh%volume_element,
m_zero, real64) , &
244 ubound(hpsi, dim = 1), &
246 ubound(psi, dim = 1), &
249 ubound(hmss, dim = 1))
251 hmss(1:nst,1:nst) = conjg(hmss(1:nst,1:nst))
254 do in = -forder, forder
255 do im = -forder, forder
256 ii = (in + forder) * nst
257 jj = (im + forder) * nst
258 hfloquet(ik_count, ii+1:ii+nst, jj+1:jj+nst) = &
259 hfloquet(ik_count, ii+1:ii+nst, jj+1:jj+nst) &
260 + hmss(1:nst, 1:nst) *
exp(-(in - im)*
m_zi * omega * it * dt)
264 hfloquet(ik_count, ii+ist, ii+ist) = hfloquet(ik_count, ii+ist, ii+ist) + in * omega
272 hfloquet(:,:,:) =
m_one / nt * hfloquet(:,:,:)
275 if (downfolding)
then
277 safe_allocate(hfloq_eff(1:nst,1:nst))
278 safe_allocate(eigenval(1:nst))
279 safe_allocate(bands(1:nik,1:nst))
281 hfloq_eff(1:nst,1:nst) =
m_zero
287 hfloq_eff(1:nst,1:nst) = hfloquet(ik,n0+1:n0+nst,m0+1:m0+nst) + &
288 m_one/omega*(matmul(hfloquet(ik,1:nst,m0+1:m0+nst), hfloquet(ik,n1+1:n1+nst,m0+1:m0+nst))- &
289 matmul(hfloquet(ik,n1+1:n1+nst,m0+1:m0+nst), hfloquet(ik,1:nst,m0+1:m0+nst)))
293 bands(ik,1:nst) = eigenval(1:nst)
295 safe_deallocate_a(hfloq_eff)
298 safe_allocate(eigenval(1:nst*fdim))
299 safe_allocate(bands(1:nik,1:nst*fdim))
300 safe_allocate(temp(1:nst*fdim, 1:nst*fdim))
303 temp(1:nst*fdim,1:nst*fdim) = hfloquet(ik,1:nst*fdim,1:nst*fdim)
305 bands(ik,1:nst*fdim) = eigenval(1:nst*fdim)
310 if (downfolding)
then
312 filename=
"downfolded_floquet_bands"
315 filename=
"floquet_bands"
320 file =
io_open(filename, sys%namespace, action =
'write')
323 write(file,
'(e13.6, 1x)',advance=
'no') bands(ik,ist)
330 if (.not. downfolding)
then
333 bands(1:nik,1:nst*fdim) =
m_zero
335 temp(1:nst*fdim, 1:nst*fdim) =
m_zero
338 temp(ii+1:ii+nst, ii+1:ii+nst) = hfloquet(ik, ii+1:ii+nst, ii+1:ii+nst)
341 bands(ik, 1:nst*fdim) = eigenval(1:nst*fdim)
345 filename=
'trivial_floquet_bands'
346 file =
io_open(filename, sys%namespace, action =
'write')
349 write(file,
'(e13.6, 1x)',advance=
'no') bands(ik,ist)
360 safe_deallocate_a(hmss)
361 safe_deallocate_a(psi)
362 safe_deallocate_a(hpsi)
363 safe_deallocate_a(temp_state1)
364 safe_deallocate_a(hfloquet)
365 safe_deallocate_a(eigenval)
366 safe_deallocate_a(bands)
367 safe_deallocate_a(temp)
subroutine floquet_init()
subroutine floquet_solve_non_interacting()
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
double exp(double __x) __attribute__((__nothrow__
This module contains interfaces for BLAS routines You should not use these routines directly....
This module handles the calculation mode.
type(calc_mode_par_t), public calc_mode_par
Singleton instance of parallel calculation mode.
integer, parameter, public p_strategy_states
parallelization in states
This module implements a calculator for the density and defines related functions.
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
subroutine, public fft_all_init(namespace)
initialize the table
subroutine, public fft_all_end()
delete all plans
real(real64), parameter, public m_two
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
real(real64), parameter, public m_zero
real(real64), parameter, public m_pi
some mathematical constants
complex(real64), parameter, public m_z0
complex(real64), parameter, public m_zi
real(real64), parameter, public m_epsilon
subroutine, public global_init(communicator)
Initialise Octopus.
real(real64), parameter, public m_one
This module implements the underlying real-space grid.
subroutine, public hamiltonian_elec_epot_generate(this, namespace, space, gr, ions, ext_partners, st, time)
subroutine, public zhamiltonian_elec_apply_all(hm, namespace, gr, st, hst)
subroutine, public io_init(defaults)
If the argument defaults is present and set to true, then the routine will not try to read anything f...
subroutine, public io_close(iunit, grp)
subroutine, public io_end()
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
This module defines the meshes, which are used in Octopus.
subroutine, public messages_end()
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
character(len=512), private msg
subroutine, public messages_init(output_dir)
subroutine, public print_date(str)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_experimental(name, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
type(mpi_grp_t), public mpi_world
This module handles the communicators for the various parallelization strategies.
type(namespace_t), public global_namespace
subroutine, public parser_init()
Initialise the Octopus parser.
subroutine, public parser_end()
End the Octopus parser.
subroutine, public profiling_end(namespace)
subroutine, public profiling_init(namespace)
Create profiling subdirectory.
integer, parameter, public restart_gs
integer, parameter, public restart_type_load
logical function, public state_kpt_is_local(st, ist, ik)
check whether a given state (ist, ik) is on the local node
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
This module handles reading and writing restart information for the states_elec_t.
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...
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
This module defines the unit system, used for input and output.
subroutine, public unit_system_init(namespace)
type(unit_system_t), public units_inp
the units systems for reading and writing
This module is intended to contain simple general-purpose utility functions and procedures.
subroutine, public print_header()
This subroutine prints the logo followed by information about the compilation and the system....
subroutine, public v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval, time, calc_energy, calc_current, force_semilocal)
Class describing the electron system.
The states_elec_t class contains all electronic wave functions.