Octopus
target_tdlocal.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
20#include "global.h"
21
24 use debug_oct_m
26 use global_oct_m
27 use grid_oct_m
28 use io_oct_m
29 use ions_oct_m
32 use, intrinsic :: iso_fortran_env
33 use mesh_oct_m
38 use parser_oct_m
40 use space_oct_m
42 use string_oct_m
44 use td_oct_m
45 use unit_oct_m
47
48
49 implicit none
50
51 private
52 public :: &
60
61
62contains
63
64 ! ----------------------------------------------------------------------
66 subroutine target_init_tdlocal(gr, namespace, tg, td)
67 type(grid_t), intent(in) :: gr
68 type(namespace_t), intent(in) :: namespace
69 type(target_t), intent(inout) :: tg
70 type(td_t), intent(in) :: td
71
72 type(block_t) :: blk
73 push_sub(target_init_tdlocal)
74
75 tg%move_ions = ion_dynamics_ions_move(td%ions_dyn)
76 tg%dt = td%dt
77
78 !%Variable OCTTdTarget
79 !%Type block
80 !%Section Calculation Modes::Optimal Control
81 !%Description
82 !% (Experimental) If <tt>OCTTargetOperator = oct_tg_td_local</tt>, then you must supply
83 !% a OCTTdTarget block. The block should only contain one element, a string cotaining the
84 !% definition of the time-dependent local target, <i>i.e.</i> a function of x,y,z and t that
85 !% is to be maximized along the evolution.
86 !%End
87 if (parse_block(namespace, 'OCTTdTarget', blk) == 0) then
88 call parse_block_string(blk, 0, 0, tg%td_local_target)
89 call conv_to_c_string(tg%td_local_target)
90 safe_allocate(tg%rho(1:gr%np))
91 call parse_block_end(blk)
92 else
93 message(1) = 'If OCTTargetOperator = oct_tg_td_local, you must supply a OCTTdTarget block.'
94 call messages_fatal(1, namespace=namespace)
95 end if
96 safe_allocate(tg%td_fitness(0:td%max_iter))
97 tg%td_fitness = m_zero
98 call target_build_tdlocal(tg, gr, m_zero)
99
100 pop_sub(target_init_tdlocal)
101 end subroutine target_init_tdlocal
102
103
104 ! ----------------------------------------------------------------------
106 subroutine target_end_tdlocal(tg)
107 type(target_t), intent(inout) :: tg
108
109 push_sub(target_end_tdlocal)
110
111 safe_deallocate_a(tg%rho)
112 safe_deallocate_a(tg%td_fitness)
113
114 pop_sub(target_end_tdlocal)
115 end subroutine target_end_tdlocal
116
117
118 ! ----------------------------------------------------------------------
119 subroutine target_output_tdlocal(tg, namespace, space, gr, dir, ions, outp)
120 type(target_t), intent(inout) :: tg
121 type(namespace_t), intent(in) :: namespace
122 class(space_t), intent(in) :: space
123 type(grid_t), intent(in) :: gr
124 character(len=*), intent(in) :: dir
125 type(ions_t), intent(in) :: ions
126 type(output_t), intent(in) :: outp
127
128 integer :: ierr
129 push_sub(target_output_tdlocal)
130
131 call io_mkdir(trim(dir), namespace)
132 call target_build_tdlocal(tg, gr, m_zero)
133 call dio_function_output(outp%how(0), trim(dir), 'td_local_target', namespace, space, gr, &
134 tg%rho, units_out%length**(-space%dim), ierr, pos=ions%pos, atoms=ions%atom)
135
136 pop_sub(target_output_tdlocal)
137 end subroutine target_output_tdlocal
138 ! ----------------------------------------------------------------------
139
140
141 ! ----------------------------------------------------------------------
143 real(real64) function target_j1_tdlocal(tg) result(j1)
144 type(target_t), intent(in) :: tg
145
146 integer :: maxiter
147 push_sub(target_j1_tdlocal)
148
149 maxiter = size(tg%td_fitness) - 1
150 j1 = m_half * tg%dt * tg%td_fitness(0) + &
151 m_half * tg%dt * tg%td_fitness(maxiter) + &
152 tg%dt * sum(tg%td_fitness(1:maxiter-1))
153
154
155 pop_sub(target_j1_tdlocal)
156 end function target_j1_tdlocal
157
158
159 ! ----------------------------------------------------------------------
161 subroutine target_chi_tdlocal(chi_out)
162 type(states_elec_t), intent(inout) :: chi_out
163
164 integer :: ik, ib
165 push_sub(target_chi_tdlocal)
166
167 !We assume that there is no time-independent operator.
168
169 do ik = chi_out%d%kpt%start, chi_out%d%kpt%end
170 do ib = chi_out%group%block_start, chi_out%group%block_end
171 call batch_set_zero(chi_out%group%psib(ib, ik))
172 end do
173 end do
174
175 pop_sub(target_chi_tdlocal)
176 end subroutine target_chi_tdlocal
177
178
179 ! ---------------------------------------------------------
182 subroutine target_tdcalc_tdlocal(tg, gr, psi, time)
183 type(target_t), intent(inout) :: tg
184 type(grid_t), intent(in) :: gr
185 type(states_elec_t), intent(in) :: psi
186 integer, intent(in) :: time
187
188 complex(real64), allocatable :: opsi(:, :), zpsi(:, :)
189 integer :: ist, ip
190 push_sub(target_tdcalc_tdlocal)
191
192 tg%td_fitness(time) = m_zero
193
194 safe_allocate(zpsi(1:gr%np, 1:psi%d%dim))
195
196 !!!! WARNING Here one should build the time-dependent target.
197 select case (psi%d%ispin)
198 case (unpolarized)
199 assert(psi%nik == 1)
200 safe_allocate(opsi(1:gr%np_part, 1))
201 opsi = m_z0
202 do ist = psi%st_start, psi%st_end
203
204 call states_elec_get_state(psi, gr, ist, 1, zpsi)
205
206 do ip = 1, gr%np
207 opsi(ip, 1) = tg%rho(ip)*zpsi(ip, 1)
208 end do
209
210 tg%td_fitness(time) = tg%td_fitness(time) + psi%occ(ist, 1)*real(zmf_dotp(gr, psi%d%dim, zpsi, opsi), real64)
211
212 end do
213 safe_deallocate_a(opsi)
214 case (spin_polarized)
215 message(1) = 'Error in target.target_tdcalc: spin_polarized.'
216 call messages_fatal(1)
217 case (spinors)
218 message(1) = 'Error in target.target_tdcalc: spinors.'
219 call messages_fatal(1)
220 end select
221
222 safe_deallocate_a(zpsi)
223
224 pop_sub(target_tdcalc_tdlocal)
225 end subroutine target_tdcalc_tdlocal
226 ! ----------------------------------------------------------------------
227
228
229 !----------------------------------------------------------
230 subroutine target_build_tdlocal(tg, gr, time)
231 type(target_t), intent(inout) :: tg
232 type(grid_t), intent(in) :: gr
233 real(real64), intent(in) :: time
234
235 integer :: ip
236 real(real64) :: xx(gr%box%dim), rr, re, im
237
238 push_sub(target_build_tdlocal)
239
240 do ip = 1, gr%np
241 call mesh_r(gr, ip, rr, coords = xx)
242 call parse_expression(re, im, gr%box%dim, xx, rr, time, tg%td_local_target)
243 tg%rho(ip) = re
244 end do
245
246 pop_sub(target_build_tdlocal)
247 end subroutine target_build_tdlocal
248 !----------------------------------------------------------
249
250end module target_tdlocal_oct_m
251
252!! Local Variables:
253!! mode: f90
254!! coding: utf-8
255!! End:
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:116
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_half
Definition: global.F90:194
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public dio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
Top-level IO routine for functions defined on the mesh.
Definition: io.F90:114
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:311
logical pure function, public ion_dynamics_ions_move(this)
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
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:414
this module contains the low-level part of the output system
Definition: output_low.F90:115
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public conv_to_c_string(str)
converts to c string
Definition: string.F90:252
subroutine, public target_output_tdlocal(tg, namespace, space, gr, dir, ions, outp)
subroutine, public target_init_tdlocal(gr, namespace, tg, td)
subroutine, public target_end_tdlocal(tg)
subroutine, public target_tdcalc_tdlocal(tg, gr, psi, time)
subroutine, public target_chi_tdlocal(chi_out)
real(real64) function, public target_j1_tdlocal(tg)
subroutine, public target_build_tdlocal(tg, gr, time)
Definition: td.F90:114
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out