Octopus
target_classical.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#include "global.h"
19
22 use debug_oct_m
23 use global_oct_m
24 use ions_oct_m
26 use, intrinsic :: iso_fortran_env
31 use parser_oct_m
34 use string_oct_m
36 use td_oct_m
37
38 implicit none
39
40 private
41 public :: &
47
48
49contains
50
51 ! ----------------------------------------------------------------------
53 subroutine target_init_classical(ions, namespace, tg, td, oct)
54 type(ions_t), intent(in) :: ions
55 type(namespace_t), intent(in) :: namespace
56 type(target_t), intent(inout) :: tg
57 type(td_t), intent(in) :: td
58 type(oct_t), intent(in) :: oct
59
60 integer :: jj, ist, jst
61 type(block_t) :: blk
62 character(len=1024) :: expression
63 push_sub(target_init_classical)
64
65 tg%move_ions = td%ions_dyn%ions_move()
66 tg%dt = td%dt
67 assert(tg%move_ions)
68
69 !%Variable OCTClassicalTarget
70 !%Type block
71 !%Section Calculation Modes::Optimal Control
72 !%Description
73 !% If <tt>OCTTargetOperator = oct_tg_classical</tt>, the you must supply this block.
74 !% It should contain a string (e.g. "(q[1,1]-q[1,2])*p[2,1]") with a mathematical
75 !% expression in terms of two arrays, q and p, that represent the position and momenta
76 !% of the classical variables. The first index runs through the various classical particles,
77 !% and the second index runs through the spatial dimensions.
78 !%
79 !% In principle, the block only contains one entry (string). However, if the expression is very
80 !% long, you can split it into various lines (one column each) that will be concatenated.
81 !%
82 !% The QOCT algorithm will attempt to maximize this expression, at the end of the propagation.
83 !%End
84 if (parse_block(namespace, 'OCTClassicalTarget', blk) == 0) then
85 tg%classical_input_string = " "
86 do jj = 0, parse_block_n(blk) - 1
87 call parse_block_string(blk, jj, 0, expression)
88 tg%classical_input_string = trim(tg%classical_input_string) // trim(expression)
89 end do
90 call parse_block_end(blk)
91 else
92 message(1) = 'If OCTTargetOperator = oct_tg_classical, then you must give the shape'
93 message(2) = 'of this target in the block "OCTClassicalTarget".'
94 call messages_fatal(2, namespace=namespace)
95 end if
96
97 !%Variable OCTMomentumDerivatives
98 !%Type block
99 !%Section Calculation Modes::Optimal Control
100 !%Description
101 !% This block should contain the derivatives of the expression given in
102 !% <tt>OCTClassicalTarget</tt> with respect to the p array components.
103 !% Each line corresponds to a different classical particle, whereas the
104 !% columns correspond to each spatial dimension: the (i,j) block component
105 !% corresponds with the derivative wrt p[i,j].
106 !%End
107 if (parse_block(namespace, 'OCTMomentumDerivatives', blk) == 0) then
108 safe_allocate(tg%mom_der_array(1:ions%natoms,1:ions%space%dim))
109 do ist = 0, ions%natoms - 1
110 do jst = 0, ions%space%dim - 1
111 call parse_block_string(blk, ist, jst, tg%mom_der_array(ist+1, jst+1))
112 end do
113 end do
114 call parse_block_end(blk)
115 else if (oct%algorithm == option__octscheme__oct_cg .or. oct%algorithm == option__octscheme__oct_bfgs) then
116 message(1) = 'If "OCTTargetOperator = oct_classical" and "OCTScheme = oct_cg" or'
117 message(2) = '"OCTScheme = oct_bfgs", then you must define the blocks "OCTClassicalTarget",'
118 message(3) = '"OCTPositionDerivatives" AND "OCTMomentumDerivatives"'
119 call messages_fatal(3, namespace=namespace)
120 end if
121
122 !%Variable OCTPositionDerivatives
123 !%Type block
124 !%Section Calculation Modes::Optimal Control
125 !%Description
126 !% This block should contain the derivatives of the expression given in
127 !% <tt>OCTClassicalTarget</tt> with respect to the q array components.
128 !% Each line corresponds to a different classical particle, whereas the
129 !% columns correspond to each spatial dimension: the (i,j) block component
130 !% corresponds with the derivative wrt q[i,j].
131 !%End
132 if (parse_block(namespace, 'OCTPositionDerivatives', blk) == 0) then
133 safe_allocate(tg%pos_der_array(1:ions%natoms,1:ions%space%dim))
134 do ist = 0, ions%natoms-1
135 do jst = 0, ions%space%dim-1
136 call parse_block_string(blk, ist, jst, tg%pos_der_array(ist+1, jst+1))
137 end do
138 end do
139 call parse_block_end(blk)
140 else if (oct%algorithm == option__octscheme__oct_cg .or. oct%algorithm == option__octscheme__oct_bfgs) then
141 message(1) = 'If "OCTTargetOperator = oct_tg_classical" and "OCTScheme = oct_cg" or'
142 message(2) = '"OCTScheme = oct_bfgs", then you must define the blocks "OCTClassicalTarget",'
143 message(3) = '"OCTPositionDerivatives" AND "OCTMomentumDerivatives"'
144 call messages_fatal(3, namespace=namespace)
145 end if
147 pop_sub(target_init_classical)
148 end subroutine target_init_classical
149 ! ----------------------------------------------------------------------
150
151 ! ----------------------------------------------------------------------
153 subroutine target_end_classical(tg)
154 type(target_t), intent(inout) :: tg
155
156 push_sub(target_end_classical)
157
158 safe_deallocate_a(tg%pos_der_array)
159 safe_deallocate_a(tg%mom_der_array)
160
161 pop_sub(target_end_classical)
162 end subroutine target_end_classical
163 ! ----------------------------------------------------------------------
164
165
166 ! ----------------------------------------------------------------------
167 subroutine target_output_classical
170 end subroutine target_output_classical
171 ! ----------------------------------------------------------------------
172
173
174 ! ----------------------------------------------------------------------
176 real(real64) function target_j1_classical(tg, qcpsi) result(j1)
177 type(target_t), intent(inout) :: tg
178 type(opt_control_state_t), intent(in) :: qcpsi
179
180 real(real64), pointer :: q(:, :), p(:, :)
181 real(real64) :: dummy(3)
182 character(len=4096) :: inp_string
183
184 push_sub(target_j1_classical)
185
186 q => opt_control_point_q(qcpsi)
187 p => opt_control_point_p(qcpsi)
188
189 inp_string = tg%classical_input_string
190 call parse_array(inp_string, q, 'q')
191 call parse_array(inp_string, p, 'p')
192 call conv_to_c_string(inp_string)
193 call parse_expression(j1, dummy(1), 1, dummy(1:3), dummy(1), dummy(1), inp_string)
194
195 nullify(q)
196 pop_sub(target_j1_classical)
197 end function target_j1_classical
198 ! ----------------------------------------------------------------------
199
200
201 ! ----------------------------------------------------------------------
203 subroutine target_chi_classical(tg, qcpsi, qcchi, ions)
204 type(target_t), intent(inout) :: tg
205 type(opt_control_state_t), intent(inout) :: qcpsi
206 type(opt_control_state_t), intent(inout) :: qcchi
207 type(ions_t), intent(in) :: ions
208
209 integer :: ist, jst, ib, iqn
210 character(len=1024) :: temp_string
211 real(real64) :: df_dv, dummy(3)
212 real(real64), pointer :: q(:, :), p(:, :), tq(:, :), tp(:, :)
213 type(states_elec_t), pointer :: chi
214 push_sub(target_chi_classical)
215
216 tq => opt_control_point_q(qcchi)
217 tp => opt_control_point_p(qcchi)
218 q => opt_control_point_q(qcpsi)
219 p => opt_control_point_p(qcpsi)
220 tq = m_zero
221 tp = m_zero
222
223 do ist = 1, ions%natoms
224 do jst=1, ions%space%dim
225 temp_string = tg%mom_der_array(ist, jst)
226 call parse_array(temp_string, p, 'p')
227 call parse_array(temp_string, q, 'q')
228 call conv_to_c_string(temp_string)
229 call parse_expression(df_dv, dummy(1), 1, dummy(1:3), dummy(1), dummy(1), temp_string)
230 tq(ist, jst) = -df_dv
231 end do
232 end do
233
234 do ist = 1, ions%natoms
235 do jst=1, ions%space%dim
236 temp_string = tg%pos_der_array(ist, jst)
237 call parse_array(temp_string, p, 'p')
238 call parse_array(temp_string, q, 'q')
239 call conv_to_c_string(temp_string)
240 call parse_expression(df_dv, dummy(1), 1, dummy(1:3), dummy(1), dummy(1), temp_string)
241 tp(ist, jst) = df_dv
242 end do
243 end do
244
245 chi => opt_control_point_qs(qcchi)
247 !We assume that there is no time-independent operator.
248 do iqn = chi%d%kpt%start, chi%d%kpt%end
249 do ib = chi%group%block_start, chi%group%block_end
250 call batch_set_zero(chi%group%psib(ib, iqn))
251 end do
252 end do
253
254 nullify(tq)
255 nullify(tp)
256 nullify(q)
257 nullify(p)
258 nullify(chi)
259 pop_sub(target_chi_classical)
260 end subroutine target_chi_classical
261 ! ----------------------------------------------------------------------
262
263end module target_classical_oct_m
264!! Local Variables:
265!! mode: f90
266!! coding: utf-8
267!! End:
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:116
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 definition of the oct_t data type, which contains some of the basic informat...
This module holds the "opt_control_state_t" datatype, which contains a quantum-classical state.
real(real64) function, dimension(:, :), pointer, public opt_control_point_p(ocs)
real(real64) function, dimension(:, :), pointer, public opt_control_point_q(ocs)
subroutine, public parse_array(inp_string, x, arraychar)
A very primitive way to "preprocess" a string that contains reference to the elements of a two-dimens...
Definition: parser.F90:700
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
real(real64) function, public target_j1_classical(tg, qcpsi)
subroutine, public target_chi_classical(tg, qcpsi, qcchi, ions)
subroutine, public target_init_classical(ions, namespace, tg, td, oct)
subroutine, public target_end_classical(tg)
subroutine, public target_output_classical
Definition: td.F90:114
This is the datatype that contains the objects that are propagated: in principle this could be both t...