Octopus
target.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 target_oct_m
22 use batch_oct_m
24 use debug_oct_m
27 use epot_oct_m
30 use fft_oct_m
31 use global_oct_m
32 use grid_oct_m
35 use io_oct_m
38 use ions_oct_m
39 use, intrinsic :: iso_fortran_env
43 use lasers_oct_m
44 use loct_oct_m
45 use math_oct_m
46 use mesh_oct_m
53 use output_oct_m
55 use parser_oct_m
58 use space_oct_m
65 use string_oct_m
67 use td_oct_m
68 use types_oct_m
69 use unit_oct_m
72
73 implicit none
74
75 private
76 public :: &
77 target_t, &
80 target_end, &
83 target_inh, &
86 target_j1, &
87 target_chi, &
91
92
93 integer, public, parameter :: &
94 oct_tg_groundstate = 1, &
95 oct_tg_excited = 2, &
98 oct_tg_jdensity = 5, &
99 oct_tg_local = 6, &
100 oct_tg_td_local = 7, &
102 oct_tg_hhg = 9, &
103 oct_tg_velocity = 10, &
104 oct_tg_hhgnew = 12, &
105 oct_tg_classical = 13, &
106 oct_tg_spin = 14
107
108 integer, public, parameter :: &
109 oct_targetmode_static = 0, &
111
112 integer, public, parameter :: &
113 oct_no_curr = 0, &
115 oct_max_curr_ring = 2, &
117
118 type target_t
119 private
120 integer :: type
121 type(states_elec_t) :: st
122 type(excited_states_t) :: est
123 real(real64), allocatable :: rho(:)
124 real(real64), allocatable :: td_fitness(:)
125 character(len=200) :: td_local_target
126 character(len=80) :: excluded_states_list
127 character(len=4096) :: vel_input_string
128 character(len=4096) :: classical_input_string
129 character(len=1024), allocatable :: vel_der_array(:,:)
130 character(len=1024), allocatable :: mom_der_array(:,:)
131 character(len=1024), allocatable :: pos_der_array(:,:)
132 real(real64), allocatable :: grad_local_pot(:,:,:)
133 logical :: move_ions
134 integer :: hhg_nks
135 integer, allocatable :: hhg_k(:)
136 real(real64), allocatable :: hhg_alpha(:)
137 real(real64), allocatable :: hhg_a(:)
138 real(real64) :: hhg_w0
139 real(real64) :: dt
140 integer :: curr_functional
141 real(real64) :: density_weight
142 real(real64) :: curr_weight
143 integer :: strt_iter_curr_tg
144 real(real64), allocatable :: spatial_curr_wgt(:)
145 character(len=1000) :: plateau_string
146 complex(real64), allocatable :: acc(:, :)
147 complex(real64), allocatable :: vel(:, :)
148 complex(real64), allocatable :: gvec(:, :)
149 real(real64), allocatable :: alpha(:)
150 complex(real64) :: spin_matrix(2, 2)
151 type(fft_t) :: fft_handler
152 end type target_t
153
154contains
155
156 ! ---------------------------------------------------------
161 subroutine target_init_propagation(tg)
162 type(target_t), intent(inout) :: tg
164
165 select case (tg%type)
166 case (oct_tg_hhgnew)
167 tg%vel = m_z0
168 tg%gvec = m_z0
169 tg%acc = m_z0
170 end select
171
173 end subroutine target_init_propagation
174
175
176 ! ----------------------------------------------------------------------
178 subroutine target_get_state(tg, st)
179 type(target_t), intent(in) :: tg
180 type(states_elec_t), intent(inout) :: st
181
182 push_sub(target_get_state)
183 call states_elec_copy(st, tg%st)
184
185 pop_sub(target_get_state)
186 end subroutine target_get_state
187 ! ----------------------------------------------------------------------
188
189
190 ! ----------------------------------------------------------------------
192 subroutine target_init(gr, kpoints, namespace, space, ions, qcs, td, w0, tg, oct, ep, mc)
193 type(grid_t), intent(in) :: gr
194 type(kpoints_t), intent(in) :: kpoints
195 type(namespace_t), intent(in) :: namespace
196 class(space_t), intent(in) :: space
197 type(ions_t), intent(in) :: ions
198 type(opt_control_state_t), intent(inout) :: qcs
199 type(td_t), intent(in) :: td
200 real(real64), intent(in) :: w0
201 type(target_t), intent(inout) :: tg
202 type(oct_t), intent(in) :: oct
203 type(epot_t), intent(inout) :: ep
204 type(multicomm_t), intent(in) :: mc
206 integer :: ierr
207 type(states_elec_t), pointer :: stin
208 type(restart_t) :: restart
209
210 push_sub(target_init)
212 stin => opt_control_point_qs(qcs)
214 !%Variable OCTTargetOperator
215 !%Type integer
216 !%Section Calculation Modes::Optimal Control
217 !%Default oct_tg_gstransformation
218 !%Description
219 !% The variable <tt>OCTTargetOperator</tt> prescribes which kind of target functional is
220 !% to be used.
221 !%Option oct_tg_groundstate 1
222 !% The target operator is a projection operator on the ground state, <i>i.e.</i> the
223 !% objective is to populate the ground state as much as possible.
224 !%Option oct_tg_excited 2
225 !% (Experimental) The target operator is an "excited state". This means that the target operator
226 !% is a linear combination of Slater determinants, each one formed by replacing
227 !% in the ground-state Slater determinant one occupied state with one excited
228 !% state (<i>i.e.</i> "single excitations"). The description of which excitations are
229 !% used, and with which weights, should be given in a file called
230 !% <tt>oct-excited-state-target</tt>.
231 !% See the documentation of subroutine <tt>excited_states_elec_init</tt> in the source
232 !% code in order to use this feature.
233 !%Option oct_tg_gstransformation 3
234 !% The target operator is a projection operator on a transformation of the ground-state
235 !% orbitals defined by the block <tt>OCTTargetTransformStates</tt>.
236 !%Option oct_tg_userdefined 4
237 !% (Experimental) Allows to define target state by using <tt>OCTTargetUserdefined</tt>.
238 !%Option oct_tg_jdensity 5
239 !% (Experimental)
240 !%Option oct_tg_local 6
241 !% (Experimental) The target operator is a local operator.
242 !%Option oct_tg_td_local 7
243 !% (Experimental) The target operator is a time-dependent local operator.
244 !%Option oct_tg_exclude_state 8
245 !% (Experimental) Target operator is the projection onto the complement of a given state, given by the
246 !% block <tt>OCTTargetTransformStates</tt>. This means that the target operator is the unity
247 !% operator minus the projector onto that state.
248 !%Option oct_tg_hhg 9
249 !% (Experimental) The target is the optimization of the HHG yield. You must supply the <tt>OCTOptimizeHarmonicSpectrum</tt>
250 !% block, and it attempts to optimize the maximum of the spectrum around each harmonic peak. You may
251 !% use only one of the gradient-less optimization schemes.
252 !%Option oct_tg_velocity 10
253 !% (Experimental) The target is a function of the velocities of the nuclei at the end of the influence of
254 !% the external field, defined by <tt>OCTVelocityTarget</tt>
255 !%Option oct_tg_hhgnew 12
256 !% (Experimental) The target is the optimization of the HHG yield. You must supply the
257 !% <tt>OCTHarmonicWeight</tt> string. It attempts to optimize the integral of the harmonic spectrum multiplied
258 !% by some user-defined weight function.
259 !%Option oct_tg_classical 13
260 !% (Experimental)
261 !%Option oct_tg_spin 14
262 !% (Experimental)
263 !%End
264 call parse_variable(namespace, 'OCTTargetOperator', oct_tg_gstransformation, tg%type)
265 if (tg%type == oct_tg_excited) call messages_experimental('OCTTargetOperator = oct_tg_excited', namespace=namespace)
266 if (tg%type == oct_tg_userdefined) call messages_experimental('OCTTargetOperator = oct_tg_userdefined')
267 if (tg%type == oct_tg_jdensity) call messages_experimental('OCTTargetOperator = oct_tg_jdensity', namespace=namespace)
268 if (tg%type == oct_tg_local) call messages_experimental('OCTTargetOperator = oct_tg_local', namespace=namespace)
269 if (tg%type == oct_tg_td_local) call messages_experimental('OCTTargetOperator = oct_tg_td_local', namespace=namespace)
270 if (tg%type == oct_tg_exclude_state) call messages_experimental('OCTTargetOperator = oct_tg_exclude_state', namespace=namespace)
271 if (tg%type == oct_tg_hhg) call messages_experimental('OCTTargetOperator = oct_tg_hhg', namespace=namespace)
272 if (tg%type == oct_tg_velocity) call messages_experimental('OCTTargetOperator = oct_tg_velocity', namespace=namespace)
273 if (tg%type == oct_tg_hhgnew) call messages_experimental('OCTTargetOperator = oct_tg_hhgnew', namespace=namespace)
274 if (tg%type == oct_tg_classical) call messages_experimental('OCTTargetOperator = oct_tg_classical', namespace=namespace)
275 if (tg%type == oct_tg_spin) call messages_experimental('OCTTargetOperator = oct_tg_spin', namespace=namespace)
276
277
278 if (.not. varinfo_valid_option('OCTTargetOperator', tg%type)) then
279 call messages_input_error(namespace, 'OCTTargetOperator')
280 end if
281
282 call states_elec_copy(tg%st, stin)
285 call restart_init(restart, namespace, restart_gs, restart_type_load, mc, ierr, mesh=gr, exact=.true.)
286 if (ierr /= 0) then
287 message(1) = "Could not read gs for OCTTargetOperator."
288 call messages_fatal(1, namespace=namespace)
289 end if
290
291 tg%curr_functional = oct_no_curr
292
293 select case (tg%type)
294 case (oct_tg_groundstate)
295 call target_init_groundstate(gr, namespace, space, tg, td, restart, kpoints)
296 case (oct_tg_excited)
297 call messages_experimental('OCTTargetOperator = oct_tg_excited', namespace=namespace)
298 call target_init_excited(gr, namespace, space, tg, td, restart, kpoints)
300 call target_init_exclude(gr, namespace, space, tg, td, restart, kpoints)
302 call target_init_gstransformation(gr, namespace, space, tg, td, restart, kpoints)
303 case (oct_tg_userdefined)
304 call target_init_userdefined(gr, namespace, tg, td)
305 case (oct_tg_jdensity)
306 call target_init_density(gr, kpoints, namespace, space, tg, stin, td, restart)
307 case (oct_tg_local)
308 call target_init_local(gr, namespace, tg, td)
309 case (oct_tg_td_local)
310 call target_init_tdlocal(gr, namespace, tg, td)
311 case (oct_tg_hhg)
312 call target_init_hhg(tg, namespace, td, w0)
313 case (oct_tg_hhgnew)
314 call messages_experimental('OCTTargetOperator = oct_tg_hhgnew', namespace=namespace)
315 call target_init_hhgnew(gr, namespace, tg, td, ions, ep)
316 case (oct_tg_velocity)
317 call target_init_velocity(gr, namespace, ions, tg, oct, td, ep)
318 case (oct_tg_classical)
319 call messages_experimental('OCTTargetOperator = oct_tg_classical', namespace=namespace)
320 call target_init_classical(ions, namespace, tg, td, oct)
321 case (oct_tg_spin)
322 call messages_experimental('OCTTargetOperator = oct_tg_spin', namespace=namespace)
323 call target_init_spin(tg, namespace)
324 case default
325 write(message(1),'(a)') "Target Operator not properly defined."
326 call messages_fatal(1, namespace=namespace)
327 end select
328
329 call restart_end(restart)
330
331 nullify(stin)
332 pop_sub(target_init)
333 end subroutine target_init
334 ! ----------------------------------------------------------------------
335
336
337 ! ----------------------------------------------------------------------
338 subroutine target_end(tg, oct)
339 type(target_t), intent(inout) :: tg
340 type(oct_t), intent(in) :: oct
341
342 push_sub(target_end)
343
344 call states_elec_end(tg%st)
345
346 select case (tg%type)
347 case (oct_tg_groundstate)
349 case (oct_tg_excited)
350 call target_end_excited()
352 call target_end_exclude()
355 case (oct_tg_userdefined)
357 case (oct_tg_jdensity)
358 call target_end_density(tg)
359 case (oct_tg_local)
360 call target_end_local(tg)
361 case (oct_tg_td_local)
362 call target_end_tdlocal(tg)
363 case (oct_tg_hhg)
364 call target_end_hhg(tg)
365 case (oct_tg_hhgnew)
366 call target_end_hhgnew(tg, oct)
367 case (oct_tg_velocity)
368 call target_end_velocity(tg, oct)
369 case (oct_tg_classical)
370 call target_end_classical(tg)
371 end select
372
373 pop_sub(target_end)
374 end subroutine target_end
375 ! ----------------------------------------------------------------------
376
377
378 ! ----------------------------------------------------------------------
379 subroutine target_output(tg, namespace, space, gr, dir, ions, hm, outp)
380 type(target_t), intent(inout) :: tg
381 type(namespace_t), intent(in) :: namespace
382 class(space_t), intent(in) :: space
383 type(grid_t), intent(in) :: gr
384 character(len=*), intent(in) :: dir
385 type(ions_t), intent(in) :: ions
386 type(hamiltonian_elec_t), intent(in) :: hm
387 type(output_t), intent(in) :: outp
388
389 push_sub(target_output)
390
391 select case (tg%type)
392 case (oct_tg_groundstate)
393 call target_output_groundstate(tg, namespace, space, gr, dir, ions, hm, outp)
394 case (oct_tg_excited)
395 call target_output_excited(tg, namespace, space, gr, dir, ions, hm, outp)
397 call target_output_exclude(tg, namespace, space, gr, dir, ions, hm, outp)
399 call target_output_gstransformation(tg, namespace, space, gr, dir, ions, hm, outp)
400 case (oct_tg_userdefined)
401 call target_output_userdefined(tg, namespace, space, gr, dir, ions, hm, outp)
402 case (oct_tg_jdensity)
403 call target_output_density(tg, namespace, space, gr, dir, ions, outp)
404 case (oct_tg_local)
405 call target_output_local(tg, namespace, space, gr, dir, ions, outp)
406 case (oct_tg_td_local)
407 call target_output_tdlocal(tg, namespace, space, gr, dir, ions, outp)
408 case (oct_tg_hhg)
409 call target_output_hhg(tg, namespace, space, gr, dir, ions, hm, outp)
410 case (oct_tg_hhgnew)
411 call target_output_hhg(tg, namespace, space, gr, dir, ions, hm, outp)
412 case (oct_tg_velocity)
413 call target_output_velocity(tg, namespace, space, gr, dir, ions, hm, outp)
414 case (oct_tg_classical)
416 end select
417
418 pop_sub(target_output)
419 end subroutine target_output
420 ! ----------------------------------------------------------------------
421
422
423 ! ---------------------------------------------------------
427 subroutine target_tdcalc(tg, namespace, space, hm, gr, ions, ext_partners, psi, time, max_time)
428 type(target_t), intent(inout) :: tg
429 type(namespace_t), intent(in) :: namespace
430 class(space_t), intent(in) :: space
431 type(hamiltonian_elec_t), intent(inout) :: hm
432 type(grid_t), intent(in) :: gr
433 type(ions_t), intent(inout) :: ions
434 type(partner_list_t), intent(in) :: ext_partners
435 type(states_elec_t), intent(inout) :: psi
436 integer, intent(in) :: time
437 integer, intent(in) :: max_time
438
439 if (target_mode(tg) /= oct_targetmode_td) return
440
441 push_sub(target_tdcalc)
442
443 tg%td_fitness(time) = m_zero
444
445 select case (tg%type)
446 case (oct_tg_hhgnew)
447 call target_tdcalc_hhgnew(tg, gr, psi, time, max_time)
448 case (oct_tg_velocity)
449 call target_tdcalc_velocity(tg, hm, gr, ions, psi, time, max_time)
450 case (oct_tg_td_local)
451 call target_tdcalc_tdlocal(tg, gr, psi, time)
452 case (oct_tg_hhg)
453 call target_tdcalc_hhg(tg, namespace, space, hm, gr, ions, ext_partners, psi, time)
454 case (oct_tg_jdensity)
455 call target_tdcalc_density(tg, gr, hm%kpoints, psi, time)
456 case default
457 message(1) = 'Error in target.target_tdcalc: default.'
458 call messages_fatal(1, namespace=namespace)
459 end select
460
461 pop_sub(target_tdcalc)
462 end subroutine target_tdcalc
463 ! ----------------------------------------------------------------------
464
465
466
467 ! ---------------------------------------------------------------
470 subroutine target_inh(psi, gr, kpoints, tg, time, inh, iter)
471 type(states_elec_t), intent(inout) :: psi
472 type(grid_t), intent(in) :: gr
473 type(kpoints_t), intent(in) :: kpoints
474 type(target_t), intent(inout) :: tg
475 real(real64), intent(in) :: time
476 type(states_elec_t), intent(inout) :: inh
477 integer, intent(in) :: iter
478
479 integer :: ik, ist, ip, idim, ib
480 complex(real64), allocatable :: zpsi(:)
481 complex(real64) :: gvec(gr%box%dim)
482
483 push_sub(target_inh)
484
485 safe_allocate(zpsi(1:gr%np))
486
487 select case (tg%type)
488 case (oct_tg_td_local)
489
490 call target_build_tdlocal(tg, gr, time)
491
492 do ik = inh%d%kpt%start, inh%d%kpt%end
493 do ist = inh%st_start, inh%st_end
494 do idim = 1, inh%d%dim
495 call states_elec_get_state(psi, gr, idim, ist, ik, zpsi)
496 zpsi(1:gr%np) = -psi%occ(ist, ik)*tg%rho(1:gr%np)*zpsi(1:gr%np)
497 call states_elec_set_state(inh, gr, idim, ist, ik, zpsi)
498 end do
499 end do
500 end do
501
502 case (oct_tg_hhgnew)
503 gvec(:) = real(tg%gvec(iter + 1, :), real64)
504
505 do ik = inh%d%kpt%start, inh%d%kpt%end
506 do ist = inh%st_start, inh%st_end
507 do idim = 1, inh%d%dim
508 call states_elec_get_state(psi, gr, idim, ist, ik, zpsi)
509 do ip = 1, gr%np
510 zpsi(ip) = -psi%occ(ist, ik)*m_two*sum(tg%grad_local_pot(1, ip, 1:gr%box%dim)*gvec(:))*zpsi(ip)
511 end do
512 call states_elec_set_state(inh, gr, idim, ist, ik, zpsi)
513 end do
514 end do
515 end do
516
517 case (oct_tg_velocity)
518
519 do ik = inh%d%kpt%start, inh%d%kpt%end
520 do ist = inh%st_start, inh%st_end
521 do idim = 1, inh%d%dim
522 call states_elec_get_state(psi, gr, idim, ist, ik, zpsi)
523 do ip = 1, gr%np
524 zpsi(ip) = -psi%occ(ist, ik)*tg%rho(ip)*zpsi(ip)
525 end do
526 call states_elec_set_state(inh, gr, idim, ist, ik, zpsi)
527 end do
528 end do
529 end do
530
531 case (oct_tg_jdensity)
532
533 do ik = inh%d%kpt%start, inh%d%kpt%end
534 do ib = inh%group%block_start, inh%group%block_end
535 call batch_set_zero(inh%group%psib(ib, ik))
536 end do
537 end do
538
539 if (abs(nint(time/tg%dt)) >= tg%strt_iter_curr_tg) then
540 call chi_current(tg, gr, kpoints, -1.0_real64, psi, inh)
541 end if
542
543 case default
544 write(message(1),'(a)') 'Internal error in target_inh'
545 call messages_fatal(1)
546
547 end select
548
549 safe_deallocate_a(zpsi)
550
551 pop_sub(target_inh)
552 end subroutine target_inh
553 !----------------------------------------------------------
554
555
556 ! ---------------------------------------------------------
561 real(real64) function target_j1(tg, namespace, gr, kpoints, qcpsi, ions) result(j1)
562 type(target_t), intent(inout) :: tg
563 type(namespace_t), intent(in) :: namespace
564 type(grid_t), intent(in) :: gr
565 type(kpoints_t), intent(in) :: kpoints
566 type(opt_control_state_t), intent(inout) :: qcpsi
567 type(ions_t), optional, intent(in) :: ions
568
569 type(states_elec_t), pointer :: psi
570
571 psi => opt_control_point_qs(qcpsi)
572
573 push_sub(target_j1)
574
575 j1 = m_zero
576 select case (tg%type)
577 case (oct_tg_groundstate)
578 j1 = target_j1_groundstate(tg, gr, psi)
579 case (oct_tg_excited)
580 j1 = target_j1_excited(tg, namespace, gr, psi)
582 j1 = target_j1_gstransformation(tg, gr, psi)
583 case (oct_tg_userdefined)
584 j1 = target_j1_userdefined(tg, gr, psi)
585 case (oct_tg_jdensity)
586 j1 = target_j1_density(gr, kpoints, tg, psi)
587 case (oct_tg_local)
588 j1 = target_j1_local(gr, tg, psi)
589 case (oct_tg_td_local)
590 j1 = target_j1_tdlocal(tg)
592 j1 = target_j1_exclude(gr, tg, psi)
593 case (oct_tg_hhg)
594 j1 = target_j1_hhg(tg, namespace)
595 case (oct_tg_hhgnew)
596 j1 = target_j1_hhgnew(gr, tg)
597 case (oct_tg_velocity)
598 j1 = target_j1_velocity(tg, ions)
599 case (oct_tg_classical)
600 j1 = target_j1_classical(tg, qcpsi)
601 case (oct_tg_spin)
602 j1 = target_j1_spin(tg, gr, psi)
603 end select
604
605 nullify(psi)
606 pop_sub(target_j1)
607 end function target_j1
608 ! ---------------------------------------------------------
609
610
611 ! ---------------------------------------------------------
613 subroutine target_chi(tg, namespace, gr, kpoints, qcpsi_in, qcchi_out, ions)
614 type(target_t), intent(inout) :: tg
615 type(namespace_t), intent(in) :: namespace
616 type(grid_t), intent(in) :: gr
617 type(kpoints_t), intent(in) :: kpoints
618 type(opt_control_state_t), target, intent(inout) :: qcpsi_in
619 type(opt_control_state_t), target, intent(inout) :: qcchi_out
620 type(ions_t), intent(in) :: ions
621
622 real(real64), pointer :: q(:, :), p(:, :)
623 type(states_elec_t), pointer :: psi_in, chi_out
624 push_sub(target_chi)
625
626 psi_in => opt_control_point_qs(qcpsi_in)
627 chi_out => opt_control_point_qs(qcchi_out)
628
629 select case (tg%type)
630 case (oct_tg_groundstate)
631
632 call target_chi_groundstate(tg, gr, psi_in, chi_out)
633 case (oct_tg_excited)
634 call target_chi_excited(tg, namespace, gr, psi_in, chi_out)
636 call target_chi_gstransformation(tg, gr, psi_in, chi_out)
637 case (oct_tg_userdefined)
638 call target_chi_userdefined(tg, gr, psi_in, chi_out)
639 case (oct_tg_jdensity)
640 call target_chi_density(tg, gr, kpoints, psi_in, chi_out)
641 case (oct_tg_local)
642 call target_chi_local(tg, gr, psi_in, chi_out)
643 case (oct_tg_td_local)
644 call target_chi_tdlocal(chi_out)
646 call target_chi_exclude(tg, gr, psi_in, chi_out)
647 case (oct_tg_hhg)
648 call target_chi_hhg(chi_out)
649 case (oct_tg_hhgnew)
650 call target_chi_hhg(chi_out)
651 case (oct_tg_velocity)
652 call target_chi_velocity(gr, tg, chi_out, ions)
653 case (oct_tg_classical)
654 call target_chi_classical(tg, qcpsi_in, qcchi_out, ions)
655 case (oct_tg_spin)
656 call target_chi_spin(tg, gr, psi_in, chi_out)
657 end select
658
659 ! Unless the target is "classical", the co-state classical variables are zero at time t=T.
660 if (tg%type .ne. oct_tg_classical) then
661 q => opt_control_point_q(qcchi_out)
662 p => opt_control_point_p(qcchi_out)
663 q = m_zero
664 p = m_zero
665 nullify(q)
666 nullify(p)
667 end if
668
669 nullify(psi_in)
670 nullify(chi_out)
671 pop_sub(target_chi)
672 end subroutine target_chi
673
674
675 ! ----------------------------------------------------------------------
676 integer pure function target_mode(tg)
677 type(target_t), intent(in) :: tg
678
679 select case (tg%type)
682 case default
684 end select
685
686 ! allow specific current functionals to be td
687 ! Attention: yet combined with static density target,
688 ! the total target is considered td.
689 select case (tg%curr_functional)
690 case (oct_curr_square_td)
692 end select
693
694 end function target_mode
695
696
697 ! ----------------------------------------------------------------------
698 integer pure function target_type(tg)
699 type(target_t), intent(in) :: tg
700 target_type = tg%type
701 end function target_type
702 ! ----------------------------------------------------------------------
703
704
705 !-----------------------------------------------------------------------
706 integer pure function target_curr_functional(tg)
707 type(target_t), intent(in) :: tg
708 target_curr_functional = tg%curr_functional
709 end function target_curr_functional
710 !-----------------------------------------------------------------------
711
712
713 ! ----------------------------------------------------------------------
714 logical pure function target_move_ions(tg)
715 type(target_t), intent(in) :: tg
716 target_move_ions = tg%move_ions
717 end function target_move_ions
718 ! ----------------------------------------------------------------------
719
720 ! ----------------------------------------------------------------------
721 logical pure function is_spatial_curr_wgt(tg)
722 type(target_t), intent(in) :: tg
723
724 is_spatial_curr_wgt = allocated(tg%spatial_curr_wgt)
725
726 end function is_spatial_curr_wgt
727 ! ----------------------------------------------------------------------
728
729#include "target_density_inc.F90"
730#include "target_velocity_inc.F90"
731#include "target_hhg_inc.F90"
732#include "target_groundstate_inc.F90"
733#include "target_excited_inc.F90"
734#include "target_gstransformation_inc.F90"
735#include "target_exclude_inc.F90"
736#include "target_userdefined_inc.F90"
737#include "target_local_inc.F90"
738#include "target_tdlocal_inc.F90"
739#include "target_classical_inc.F90"
740#include "target_spin_inc.F90"
741
742end module target_oct_m
743
744!! Local Variables:
745!! mode: f90
746!! coding: utf-8
747!! End:
This module implements batches of mesh functions.
Definition: batch.F90:133
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:116
subroutine, public batch_set_zero(this, np, async)
fill all mesh functions of the batch with zero
Definition: batch_ops.F90:240
This module implements a calculator for the density and defines related functions.
Definition: density.F90:120
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
Definition: fft.F90:118
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_zero
Definition: global.F90:187
complex(real64), parameter, public m_z0
Definition: global.F90:197
This module implements the underlying real-space grid.
Definition: grid.F90:117
This module defines classes and functions for interaction partners.
Definition: io.F90:114
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
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:420
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1097
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
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.
type(states_elec_t) function, pointer, public opt_control_point_qs(ocs)
this module contains the output system
Definition: output_low.F90:115
this module contains the output system
Definition: output.F90:115
integer, parameter, public restart_gs
Definition: restart.F90:229
subroutine, public restart_init(restart, namespace, data_type, type, mc, ierr, mesh, dir, exact)
Initializes a restart object.
Definition: restart.F90:514
integer, parameter, public restart_type_load
Definition: restart.F90:225
subroutine, public restart_end(restart)
Definition: restart.F90:720
subroutine, public states_elec_end(st)
finalize the states_elec_t object
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.
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 target_end_groundstate()
Definition: target.F90:2330
subroutine target_end_gstransformation()
Definition: target.F90:2767
real(real64) function target_j1_groundstate(tg, gr, psi)
Definition: target.F90:2360
subroutine target_chi_excited(tg, namespace, gr, psi_in, chi_out)
Definition: target.F90:2561
logical pure function is_spatial_curr_wgt(tg)
Definition: target.F90:815
subroutine target_build_tdlocal(tg, gr, time)
Definition: target.F90:3602
integer, parameter, public oct_tg_jdensity
Definition: target.F90:186
subroutine target_output_exclude(tg, namespace, space, gr, dir, ions, hm, outp)
Definition: target.F90:2943
subroutine target_end_velocity(tg, oct)
Definition: target.F90:1681
subroutine target_tdcalc_tdlocal(tg, gr, psi, time)
Definition: target.F90:3554
subroutine target_init_hhgnew(gr, namespace, tg, td, ions, ep)
Definition: target.F90:1950
subroutine target_end_hhg(tg)
Definition: target.F90:2046
subroutine target_chi_hhg(chi_out)
Definition: target.F90:2160
subroutine target_chi_classical(tg, qcpsi, qcchi, ions)
Definition: target.F90:3800
integer, parameter, public oct_targetmode_static
Definition: target.F90:201
real(real64) function target_j1_classical(tg, qcpsi)
Definition: target.F90:3773
logical pure function, public target_move_ions(tg)
Definition: target.F90:808
subroutine target_output_gstransformation(tg, namespace, space, gr, dir, ions, hm, outp)
Definition: target.F90:2776
integer, parameter, public oct_tg_userdefined
Definition: target.F90:186
subroutine target_init_velocity(gr, namespace, ions, tg, oct, td, ep)
Definition: target.F90:1539
subroutine target_chi_tdlocal(chi_out)
Definition: target.F90:3533
subroutine target_init_hhg(tg, namespace, td, w0)
Definition: target.F90:1871
integer, parameter, public oct_tg_hhgnew
Definition: target.F90:186
integer, parameter, public oct_tg_spin
Definition: target.F90:186
integer, parameter, public oct_targetmode_td
Definition: target.F90:201
subroutine, public target_init_propagation(tg)
This routine performs all the things that must be initialized prior to a forward evolution,...
Definition: target.F90:255
subroutine target_tdcalc_hhgnew(tg, gr, psi, time, max_time)
Definition: target.F90:2180
subroutine target_init_classical(ions, namespace, tg, td, oct)
Definition: target.F90:3650
integer, parameter, public oct_max_curr_ring
Definition: target.F90:205
real(real64) function, public target_j1(tg, namespace, gr, kpoints, qcpsi, ions)
Calculates the J1 functional, i.e.: in the time-independent case, or else in the time-dependent cas...
Definition: target.F90:655
subroutine target_chi_velocity(gr, tg, chi_out, ions)
Definition: target.F90:1751
subroutine target_chi_userdefined(tg, gr, psi_in, chi_out)
Definition: target.F90:3215
subroutine target_end_excited()
Definition: target.F90:2514
subroutine target_end_density(tg)
Definition: target.F90:1146
real(real64) function target_j1_spin(tg, gr, psi)
Definition: target.F90:3942
subroutine target_output_velocity(tg, namespace, space, gr, dir, ions, hm, outp)
Definition: target.F90:1699
real(real64) function target_j1_velocity(tg, ions)
Definition: target.F90:1721
subroutine target_init_density(gr, kpoints, namespace, space, tg, stin, td, restart)
Definition: target.F90:844
real(real64) function target_j1_exclude(gr, tg, psi)
Definition: target.F90:2965
subroutine target_end_local(tg)
Definition: target.F90:3324
subroutine target_init_exclude(mesh, namespace, space, tg, td, restart, kpoints)
Definition: target.F90:2897
subroutine, public target_init(gr, kpoints, namespace, space, ions, qcs, td, w0, tg, oct, ep, mc)
The target is initialized, mainly by reading from the inp file.
Definition: target.F90:286
integer, parameter, public oct_curr_square_td
Definition: target.F90:205
subroutine target_init_excited(mesh, namespace, space, tg, td, restart, kpoints)
Definition: target.F90:2459
subroutine target_chi_density(tg, gr, kpoints, psi_in, chi_out)
Definition: target.F90:1239
real(real64) function target_j1_excited(tg, namespace, gr, psi)
Definition: target.F90:2545
subroutine target_output_classical
Definition: target.F90:3764
integer, parameter, public oct_curr_square
Definition: target.F90:205
real(real64) function target_j1_density(gr, kpoints, tg, psi)
Definition: target.F90:1190
integer, parameter, public oct_tg_velocity
Definition: target.F90:186
subroutine target_output_density(tg, namespace, space, mesh, dir, ions, outp)
Definition: target.F90:1163
subroutine target_end_userdefined()
Definition: target.F90:3153
real(real64) function target_j1_userdefined(tg, gr, psi)
Definition: target.F90:3182
subroutine target_output_tdlocal(tg, namespace, space, gr, dir, ions, outp)
Definition: target.F90:3491
integer, parameter, public oct_tg_excited
Definition: target.F90:186
integer, parameter, public oct_tg_gstransformation
Definition: target.F90:186
subroutine target_chi_spin(tg, gr, psi_in, chi_out)
Definition: target.F90:3972
real(real64) function target_j1_hhgnew(gr, tg)
Definition: target.F90:2138
real(real64) function target_j1_local(mesh, tg, psi)
Definition: target.F90:3360
subroutine target_init_groundstate(mesh, namespace, space, tg, td, restart, kpoints)
Definition: target.F90:2299
integer, parameter, public oct_tg_groundstate
Definition: target.F90:186
integer, parameter, public oct_tg_td_local
Definition: target.F90:186
subroutine, public target_get_state(tg, st)
This just copies the states_elec_t variable present in target, into st.
Definition: target.F90:272
subroutine target_end_tdlocal(tg)
Definition: target.F90:3478
integer, parameter, public oct_tg_local
Definition: target.F90:186
real(real64) function target_j1_gstransformation(tg, gr, psi)
Definition: target.F90:2798
subroutine target_tdcalc_density(tg, gr, kpoints, psi, time)
Definition: target.F90:1318
integer pure function, public target_type(tg)
Definition: target.F90:792
real(real64) function target_j1_tdlocal(tg)
Definition: target.F90:3515
subroutine, public target_inh(psi, gr, kpoints, tg, time, inh, iter)
Calculates the inhomogeneous term that appears in the equation for chi, and places it into inh.
Definition: target.F90:564
subroutine target_init_local(gr, namespace, tg, td)
Definition: target.F90:3279
subroutine, public target_output(tg, namespace, space, gr, dir, ions, hm, outp)
Definition: target.F90:473
subroutine target_end_classical(tg)
Definition: target.F90:3750
real(real64) function target_j1_hhg(tg, namespace)
Definition: target.F90:2106
subroutine target_chi_gstransformation(tg, gr, psi_in, chi_out)
Definition: target.F90:2833
subroutine target_chi_local(tg, mesh, psi_in, chi_out)
Definition: target.F90:3379
subroutine, public target_tdcalc(tg, namespace, space, hm, gr, ions, ext_partners, psi, time, max_time)
Calculates, at a given point in time marked by the integer index, the integrand of the target functio...
Definition: target.F90:521
subroutine target_output_groundstate(tg, namespace, space, gr, dir, ions, hm, outp)
Definition: target.F90:2338
subroutine target_init_userdefined(gr, namespace, tg, td)
Definition: target.F90:3063
integer, parameter, public oct_tg_hhg
Definition: target.F90:186
subroutine, public target_chi(tg, namespace, gr, kpoints, qcpsi_in, qcchi_out, ions)
Calculate .
Definition: target.F90:707
subroutine, public target_end(tg, oct)
Definition: target.F90:432
subroutine target_tdcalc_hhg(tg, namespace, space, hm, gr, ions, ext_partners, psi, time)
Definition: target.F90:2248
subroutine target_output_userdefined(tg, namespace, space, gr, dir, ions, hm, outp)
Definition: target.F90:3161
integer pure function, public target_mode(tg)
Definition: target.F90:770
subroutine target_chi_exclude(tg, gr, psi_in, chi_out)
Definition: target.F90:2997
integer pure function, public target_curr_functional(tg)
Definition: target.F90:800
subroutine target_end_hhgnew(tg, oct)
Definition: target.F90:2083
subroutine target_chi_groundstate(tg, gr, psi_in, chi_out)
Definition: target.F90:2394
subroutine target_init_tdlocal(gr, namespace, tg, td)
Definition: target.F90:3438
integer, parameter, public oct_tg_classical
Definition: target.F90:186
subroutine target_end_exclude()
Definition: target.F90:2935
subroutine target_output_excited(tg, namespace, space, gr, dir, ions, hm, outp)
Definition: target.F90:2522
integer, parameter, public oct_tg_exclude_state
Definition: target.F90:186
subroutine target_output_hhg(tg, namespace, space, gr, dir, ions, hm, outp)
Definition: target.F90:2061
subroutine chi_current(tg, gr, kpoints, factor, psi_in, chi)
Definition: target.F90:1401
subroutine target_init_spin(tg, namespace)
Definition: target.F90:3888
subroutine target_tdcalc_velocity(tg, hm, gr, ions, psi, time, max_time)
Definition: target.F90:1797
subroutine target_init_gstransformation(gr, namespace, space, tg, td, restart, kpoints)
Definition: target.F90:2723
subroutine target_output_local(tg, namespace, space, mesh, dir, ions, outp)
Definition: target.F90:3336
Definition: td.F90:114
type(type_t), public type_cmplx
Definition: types.F90:134
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.
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
Stores all communicators and groups.
Definition: multicomm.F90:206
!brief The oct_t datatype stores the basic information about how the OCT run is done.
This is the datatype that contains the objects that are propagated: in principle this could be both t...
output handler class
Definition: output_low.F90:163
The states_elec_t class contains all electronic wave functions.
int true(void)