Octopus
initst.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
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
21module initst_oct_m
22 use debug_oct_m
24 use global_oct_m
25 use grid_oct_m
27 use ions_oct_m
31 use parser_oct_m
36 use string_oct_m
38 use td_oct_m
39 use v_ks_oct_m
41 use types_oct_m
42 use xc_oct_m
43
44 implicit none
45
46 private
47 public :: initial_state_init
48
49 integer, parameter :: &
50 oct_is_groundstate = 1, &
51 oct_is_excited = 2, &
54
55
56contains
57
58
59 ! ---------------------------------------------------------
60 subroutine initial_state_init(sys, qcstate)
61 type(electrons_t), intent(inout) :: sys
62 type(opt_control_state_t), target, intent(inout) :: qcstate
63
64 integer :: ik, ib, idim, inst, inik, id, is, ip, ierr, &
65 no_states, istype, freeze_orbitals
66 type(block_t) :: blk
67 real(real64) :: xx(1:sys%space%dim), rr, psi_re, psi_im
68 type(restart_t) :: restart
69 complex(real64), allocatable :: zpsi(:, :)
70
71 type(states_elec_t), pointer :: psi
72
73 push_sub(initial_state_init)
74
75 call opt_control_state_init(qcstate, sys%st, sys%ions)
76 psi => opt_control_point_qs(qcstate)
78 call states_elec_allocate_wfns(psi, sys%gr, type_cmplx)
79
80 !%Variable OCTInitialState
81 !%Type integer
82 !%Section Calculation Modes::Optimal Control
83 !%Default oct_is_groundstate
84 !%Description
85 !% Describes the initial state of the quantum system.
86 !% Possible arguments are:
87 !%Option oct_is_groundstate 1
88 !% Start in the ground state.
89 !%Option oct_is_excited 2
90 !% Currently not in use.
91 !%Option oct_is_gstransformation 3
92 !% Start in a transformation of the ground-state orbitals, as defined in the
93 !% block <tt>OCTInitialTransformStates</tt>.
94 !%Option oct_is_userdefined 4
95 !% Start in a user-defined state.
96 !%End
97 call parse_variable(sys%namespace, 'OCTInitialState', oct_is_groundstate, istype)
98 if (.not. varinfo_valid_option('OCTInitialState', istype)) call messages_input_error(sys%namespace, 'OCTInitialState')
99
100 select case (istype)
101 case (oct_is_groundstate)
102 message(1) = 'Info: Using ground state for initial state.'
103 call messages_info(1, namespace=sys%namespace)
104 call restart_init(restart, sys%namespace, restart_gs, restart_type_load, sys%mc, ierr, mesh=sys%gr, exact=.true.)
105 if (ierr == 0) then
106 call states_elec_load(restart, sys%namespace, sys%space, psi, sys%gr, sys%kpoints, &
107 fixed_occ=.true., ierr=ierr)
108 end if
109 if (ierr /= 0) then
110 message(1) = "Unable to read wavefunctions."
111 call messages_fatal(1, namespace=sys%namespace)
112 end if
113 call restart_end(restart)
115 case (oct_is_excited)
116 message(1) = 'Using an excited state as the starting state for an '
117 message(2) = 'optimal-control run is not possible yet.'
118 message(3) = 'Try using "OCTInitialState = oct_is_transformation" instead.'
119 call messages_fatal(3, namespace=sys%namespace)
120
122 message(1) = 'Info: Using superposition of states for initial state.'
123 call messages_info(1, namespace=sys%namespace)
124
125
126 !%Variable OCTInitialTransformStates
127 !%Type block
128 !%Section Calculation Modes::Optimal Control
129 !%Description
130 !% If <tt>OCTInitialState = oct_is_gstransformation</tt>, you must specify an
131 !% <tt>OCTInitialTransformStates</tt> block, in order to specify which linear
132 !% combination of the states present in <tt>restart/gs</tt> is used to
133 !% create the initial state.
134 !%
135 !% The syntax is the same as the <tt>TransformStates</tt> block.
136 !%End
137
138 if (.not. parse_is_defined(sys%namespace, "OCTInitialTransformStates")) then
139 message(1) = 'If "OCTInitialState = oct_is_gstransformation", then you must'
140 message(2) = 'supply an "OCTInitialTransformStates" block to define the transformation.'
141 call messages_fatal(2, namespace=sys%namespace)
142 end if
143
144 call restart_init(restart, sys%namespace, restart_gs, restart_type_load, sys%mc, ierr, mesh=sys%gr, exact=.true.)
145 if (ierr /= 0) then
146 message(1) = "Could not read states for OCTInitialTransformStates."
147 call messages_fatal(1, namespace=sys%namespace)
148 end if
149
150 call states_elec_transform(psi, sys%namespace, sys%space, restart, sys%gr, sys%kpoints, prefix = "OCTInitial")
151 call restart_end(restart)
152
154 message(1) = 'Info: Building user-defined initial state.'
155 call messages_info(1, namespace=sys%namespace)
156
157 !%Variable OCTInitialUserdefined
158 !%Type block
159 !%Section Calculation Modes::Optimal Control
160 !%Description
161 !% Define an initial state. Syntax follows the one of the <tt>UserDefinedStates</tt> block.
162 !% Example:
163 !%
164 !% <tt>%OCTInitialUserdefined
165 !% <br>&nbsp;&nbsp; 1 | 1 | 1 | "exp(-r^2)*exp(-i*0.2*x)"
166 !% <br>%</tt>
167 !%
168 !%End
169 if (parse_block(sys%namespace, 'OCTInitialUserdefined', blk) == 0) then
170
171 safe_allocate(zpsi(1:sys%gr%np, 1:psi%d%dim))
172
173 no_states = parse_block_n(blk)
174 do ib = 1, no_states
175 call parse_block_integer(blk, ib - 1, 0, idim)
176 call parse_block_integer(blk, ib - 1, 1, inst)
177 call parse_block_integer(blk, ib - 1, 2, inik)
178
179 ! read formula strings and convert to C strings
180 do id = 1, psi%d%dim
181 do is = 1, psi%nst
182 do ik = 1, psi%nik
183
184 ! does the block entry match and is this node responsible?
185 if (.not. (id == idim .and. is == inst .and. ik == inik &
186 .and. psi%st_start <= is .and. psi%st_end >= is)) cycle
187
188 ! parse formula string
189 call parse_block_string( &
190 blk, ib - 1, 3, psi%user_def_states(id, is, ik))
191 ! convert to C string
192 call conv_to_c_string(psi%user_def_states(id, is, ik))
193
194 do ip = 1, sys%gr%np
195 xx = sys%gr%x(ip, :)
196 rr = norm2(xx)
197
198 ! parse user-defined expressions
199 call parse_expression(psi_re, psi_im, &
200 sys%space%dim, xx, rr, m_zero, psi%user_def_states(id, is, ik))
201 ! fill state
202 zpsi(ip, id) = cmplx(psi_re, psi_im, real64)
203 end do
204 call states_elec_set_state(psi, sys%gr, id, is, ik, zpsi(:, id))
205 end do
206 end do
207 end do
208 end do
209 call parse_block_end(blk)
210 do ik = 1, psi%nik
211 do is = psi%st_start, psi%st_end
212 call states_elec_get_state(psi, sys%gr, is, ik, zpsi)
213 call zmf_normalize(sys%gr, psi%d%dim, zpsi)
214 call states_elec_set_state(psi, sys%gr, is, ik, zpsi)
215 end do
216 end do
217 safe_deallocate_a(zpsi)
218 else
219 call messages_variable_is_block(sys%namespace, 'OCTInitialUserdefined')
220 end if
221
222 case default
223 write(message(1),'(a)') "No valid initial state defined."
224 write(message(2),'(a)') "Choosing the ground state."
225 call messages_info(2, namespace=sys%namespace)
226 end select
227
228 ! Check whether we want to freeze some of the deeper orbitals.
229 call parse_variable(sys%namespace, 'TDFreezeOrbitals', 0, freeze_orbitals)
230 if (freeze_orbitals > 0) then
231 ! In this case, we first freeze the orbitals, then calculate the Hxc potential.
232 call states_elec_freeze_orbitals(psi, sys%namespace, sys%space, sys%gr, sys%mc, sys%kpoints, &
233 freeze_orbitals, family_is_mgga(sys%ks%xc_family))
234 write(message(1),'(a,i4,a,i4,a)') 'Info: The lowest', freeze_orbitals, &
235 ' orbitals have been frozen.', psi%nst, ' will be propagated.'
236 call messages_info(1, namespace=sys%namespace)
237 call density_calc(psi, sys%gr, psi%rho)
238 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, psi, sys%ions, sys%ext_partners, calc_eigenval = .true.)
239 elseif (freeze_orbitals < 0) then
240 ! This means SAE approximation. We calculate the Hxc first, then freeze all
241 ! orbitals minus one.
242 write(message(1),'(a)') 'Info: The single-active-electron approximation will be used.'
243 call messages_info(1, namespace=sys%namespace)
244 call density_calc(psi, sys%gr, psi%rho)
245 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, psi, sys%ions, sys%ext_partners, calc_eigenval = .true.)
246 call states_elec_freeze_orbitals(psi, sys%namespace, sys%space, sys%gr, sys%mc, sys%kpoints, &
247 psi%nst - 1, family_is_mgga(sys%ks%xc_family))
248 call v_ks_freeze_hxc(sys%ks)
249 call density_calc(psi, sys%gr, psi%rho)
250 else
251 ! Normal run.
252 call density_calc(psi, sys%gr, psi%rho)
253 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, psi, sys%ions, sys%ext_partners, calc_eigenval = .true.)
254 end if
255
256 pop_sub(initial_state_init)
257 end subroutine initial_state_init
258
259end module initst_oct_m
260
261!! Local Variables:
262!! mode: f90
263!! coding: utf-8
264!! End:
This module implements a calculator for the density and defines related functions.
Definition: density.F90:120
subroutine, public states_elec_freeze_orbitals(st, namespace, space, gr, mc, kpoints, n, family_is_mgga)
Calculate partial density for frozen orbitals.
Definition: density.F90:655
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Definition: density.F90:612
real(real64), parameter, public m_zero
Definition: global.F90:188
This module implements the underlying real-space grid.
Definition: grid.F90:117
integer, parameter oct_is_excited
Definition: initst.F90:142
integer, parameter oct_is_gstransformation
Definition: initst.F90:142
subroutine, public initial_state_init(sys, qcstate)
Definition: initst.F90:154
integer, parameter oct_is_userdefined
Definition: initst.F90:142
This module defines various routines, operating on mesh functions.
subroutine, public zmf_normalize(mesh, dim, psi, norm)
Normalize a mesh function psi.
subroutine, public messages_variable_is_block(namespace, name)
Definition: messages.F90:1070
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_input_error(namespace, var, details, row, column)
Definition: messages.F90:714
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
This module holds the "opt_control_state_t" datatype, which contains a quantum-classical state.
type(states_elec_t) function, pointer, public opt_control_point_qs(ocs)
subroutine, public opt_control_state_init(ocs, qstate, ions)
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
integer, parameter, public restart_gs
Definition: restart.F90:200
subroutine, public restart_init(restart, namespace, data_type, type, mc, ierr, mesh, dir, exact)
Initializes a restart object.
Definition: restart.F90:516
integer, parameter, public restart_type_load
Definition: restart.F90:245
subroutine, public restart_end(restart)
Definition: restart.F90:722
subroutine, public states_elec_deallocate_wfns(st)
Deallocates the KS wavefunctions defined within a states_elec_t structure.
subroutine, public states_elec_allocate_wfns(st, mesh, wfs_type, skip, packed)
Allocates the KS wavefunctions defined within a states_elec_t structure.
This module handles reading and writing restart information for the states_elec_t.
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 conv_to_c_string(str)
converts to c string
Definition: string.F90:229
Definition: td.F90:114
type(type_t), parameter, public type_cmplx
Definition: types.F90:134
subroutine, public v_ks_freeze_hxc(ks)
Definition: v_ks.F90:1443
subroutine, public v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval, time, calc_energy, calc_current, force_semilocal)
Definition: v_ks.F90:744
Definition: xc.F90:114
pure logical function, public family_is_mgga(family, only_collinear)
Is the xc function part of the mGGA family.
Definition: xc.F90:582
int true(void)