Octopus
opt_control.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
26 use debug_oct_m
30 use filter_oct_m
31 use global_oct_m
32 use grid_oct_m
34 use initst_oct_m
35 use ions_oct_m
36 use iso_c_binding
37 use, intrinsic :: iso_fortran_env
38 use io_oct_m
39 use lasers_oct_m
40 use loct_oct_m
41 use math_oct_m
49 use output_oct_m
50 use parser_oct_m
59 use target_oct_m
60 use td_oct_m
61
62 implicit none
63
64 private
65 public :: &
72
73
75 type(filter_t), save :: filter
76 type(oct_t), save :: oct
77 type(oct_iterator_t), save :: iterator
78 type(target_t), save :: oct_target
79 type(opt_control_state_t), save :: initial_st
80
81
83 type(controlfunction_t), save :: par_
84 type(electrons_t), pointer :: sys_
85 type(hamiltonian_elec_t), pointer :: hm_
86 type(td_t), pointer :: td_
87 real(real64), allocatable :: x_(:)
88 integer :: index_
89
90contains
91
92 ! ---------------------------------------------------------
93 subroutine opt_control_run(system)
94 class(*), intent(inout) :: system
95
96 push_sub(opt_control_run)
97
98 select type (system)
99 class is (multisystem_basic_t)
100 message(1) = "CalculationMode = opt_control not implemented for multi-system calculations"
101 call messages_fatal(1, namespace=system%namespace)
102 type is (electrons_t)
103 call opt_control_run_legacy(system)
104 end select
105
106 pop_sub(opt_control_run)
107 end subroutine opt_control_run
108
111 subroutine opt_control_run_legacy(sys)
112 type(electrons_t), target, intent(inout) :: sys
113
114 type(td_t), target :: td
115 type(controlfunction_t) :: par, par_new, par_prev
116 logical :: stop_loop
117 real(real64) :: j1
118 type(oct_prop_t) :: prop_chi, prop_psi
119 type(states_elec_t) :: psi
120
121 type(lasers_t), pointer :: lasers
122
123 push_sub(opt_control_run_legacy)
124
125 if (sys%hm%pcm%run_pcm) then
126 call messages_not_implemented("PCM for CalculationMode /= gs or td", namespace=sys%namespace)
127 end if
128
129 if (sys%kpoints%use_symmetries) then
130 call messages_experimental("KPoints symmetries with CalculationMode = opt_control", namespace=sys%namespace)
131 end if
132
133 ! Creates a directory where the optimal control stuff will be written. The name of the directory
134 ! is stored in the preprocessor macro OCT_DIR, which should be defined in src/include/global.h
135 call io_mkdir(oct_dir, sys%namespace)
136
137 ! Initializes the time propagator. Then, it forces the propagation to be self consistent, in case
138 ! the theory level is not "independent particles".
139 call td_init(td, sys%namespace, sys%space, sys%gr, sys%ions, sys%st, sys%ks, sys%hm, sys%ext_partners, sys%outp)
140 if (sys%hm%theory_level /= independent_particles) call propagator_elec_set_scf_prop(td%tr, threshold = 1.0e-14_real64)
141
142 ! Read general information about how the OCT run will be made, from inp file. "oct_read_inp" is
143 ! in the opt_control_global_oct_m module (like the definition of the oct_t data type)
144 call oct_read_inp(oct, sys%namespace)
145
146 ! Read info about, and prepare, the control functions
147 call controlfunction_mod_init(sys%ext_partners, sys%namespace, td%dt, td%max_iter, oct%mode_fixed_fluence)
148 call controlfunction_init(par, td%dt, td%max_iter)
149 call controlfunction_set(par, sys%ext_partners)
150 ! This prints the initial control parameters, exactly as described in the inp file,
151 ! that is, without applying any envelope or filter.
152 call controlfunction_write(oct_dir//'initial_laser_inp', par, sys%namespace)
154 call controlfunction_to_h(par, sys%ext_partners)
155 call messages_print_with_emphasis(msg="TD ext. fields after processing", namespace=sys%namespace)
156 lasers => list_get_lasers(sys%ext_partners)
157 if(associated(lasers)) call laser_write_info(lasers%lasers, namespace=sys%namespace)
158 call messages_print_with_emphasis(namespace=sys%namespace)
159 call controlfunction_write(oct_dir//'initial_laser', par, sys%namespace)
160
161
162 ! Startup of the iterator data type (takes care of counting iterations, stopping, etc).
163 call oct_iterator_init(iterator, sys%namespace, par)
164
165
166 ! Initialization of the propagation_oct_m module.
167 call propagation_mod_init(td%max_iter, oct%eta, oct%delta, oct%number_checkpoints, &
168 (oct%algorithm == option__octscheme__oct_zbr98), &
169 (oct%algorithm == option__octscheme__oct_cg) .or. &
170 (oct%algorithm == option__octscheme__oct_bfgs) .or. &
171 (oct%algorithm == option__octscheme__oct_nlopt_lbfgs))
173
174 ! If filters are to be used, they also have to be initialized.
175 call filter_init(td%max_iter, sys%namespace, td%dt, filter)
176 call filter_write(filter, sys%namespace)
179 ! Figure out the starting wavefunction(s), and the target.
180 call initial_state_init(sys, initial_st)
181 call target_init(sys%gr, sys%kpoints, sys%namespace, sys%space, sys%ions, initial_st, td, &
182 controlfunction_w0(par), oct_target, oct, sys%hm%ep, sys%mc)
183
184 ! Sanity checks.
185 call check_faulty_runmodes(sys, td%tr)
187
188 ! Informative output.
189 call opt_control_get_qs(psi, initial_st)
190 call output_states(sys%outp, sys%namespace, sys%space, oct_dir//'initial', psi, sys%gr, sys%ions, sys%hm, -1)
191 call target_output(oct_target, sys%namespace, sys%space, sys%gr, oct_dir//'target', sys%ions, sys%hm, sys%outp)
192 call states_elec_end(psi)
193
194
195 ! mode switcher; here is where the real run is made.
196 select case (oct%algorithm)
197 case (option__octscheme__oct_zbr98)
198 message(1) = "Info: Starting OCT iteration using scheme: ZBR98"
199 call messages_info(1, namespace=sys%namespace)
200 call scheme_zbr98()
201 case (option__octscheme__oct_wg05)
202 message(1) = "Info: Starting OCT iteration using scheme: WG05"
203 call messages_info(1, namespace=sys%namespace)
205 case (option__octscheme__oct_zr98)
206 message(1) = "Info: Starting OCT iteration using scheme: ZR98"
207 call messages_info(1, namespace=sys%namespace)
208 call scheme_mt03()
209 case (option__octscheme__oct_mt03)
210 message(1) = "Info: Starting OCT iteration using scheme: MT03"
211 call messages_info(1, namespace=sys%namespace)
212 call scheme_mt03()
213 case (option__octscheme__oct_krotov)
214 message(1) = "Info: Starting OCT iteration using scheme: KROTOV"
215 call messages_info(1, namespace=sys%namespace)
216 call scheme_mt03()
217 case (option__octscheme__oct_straight_iteration)
218 message(1) = "Info: Starting OCT iterations using scheme: STRAIGHT ITERATION"
219 call messages_info(1, namespace=sys%namespace)
221 case (option__octscheme__oct_cg)
222 message(1) = "Info: Starting OCT iterations using scheme: CONJUGATE GRADIENTS"
223 call messages_info(1, namespace=sys%namespace)
224 call scheme_cg()
225 case (option__octscheme__oct_bfgs)
226 message(1) = "Info: Starting OCT iterations using scheme: BFGS"
227 call messages_info(1, namespace=sys%namespace)
228 call scheme_cg()
229 case (option__octscheme__oct_direct)
230 message(1) = "Info: Starting OCT iterations using scheme: DIRECT OPTIMIZATION (NELDER-MEAD)"
231 call messages_info(1, namespace=sys%namespace)
232 call scheme_direct()
233 case (option__octscheme__oct_nlopt_bobyqa)
234 message(1) = "Info: Starting OCT iterations using scheme: DIRECT OPTIMIZATION (NLOPT - BOBYQA)"
235 call messages_info(1, namespace=sys%namespace)
236 call scheme_nlopt()
237 case (option__octscheme__oct_nlopt_lbfgs)
238 message(1) = "Info: Starting OCT iterations using scheme: DIRECT OPTIMIZATION (NLOPT - LBFGS)"
239 call messages_info(1, namespace=sys%namespace)
240 call scheme_nlopt()
241 case default
242 call messages_input_error(sys%namespace, 'OCTScheme')
243 end select
244
245 ! do final test run: propagate initial state with optimal field
246 call oct_finalcheck(sys, td)
247
248 ! clean up
249 call controlfunction_end(par)
250 call oct_iterator_end(iterator, sys%namespace)
251 call filter_end(filter)
252 call td_end(td)
253 call opt_control_state_end(initial_st)
254 call target_end(oct_target, oct)
256
258
259 contains
260
261
262 ! ---------------------------------------------------------
263 subroutine scheme_straight_iteration()
265
267 call controlfunction_copy(par_new, par)
268 ctr_loop: do
269 call controlfunction_copy(par_prev, par)
270 call f_striter(sys, td, par, j1)
271 stop_loop = iteration_manager(sys%namespace, j1, par_prev, par, iterator)
272 if (clean_stop(sys%mc%master_comm) .or. stop_loop) exit ctr_loop
273 end do ctr_loop
274
275 call controlfunction_end(par_new)
276 call controlfunction_end(par_prev)
278 end subroutine scheme_straight_iteration
279 ! ---------------------------------------------------------
280
281
282 ! ---------------------------------------------------------
283 subroutine scheme_mt03()
284 type(opt_control_state_t) :: psi
286
287 call opt_control_state_null(psi)
288 call opt_control_state_copy(psi, initial_st)
289 call oct_prop_init(prop_chi, sys%namespace, "chi", sys%gr, sys%mc)
290 call oct_prop_init(prop_psi, sys%namespace, "psi", sys%gr, sys%mc)
291
292 call controlfunction_copy(par_new, par)
293 ctr_loop: do
294 call controlfunction_copy(par_prev, par)
295 call f_iter(sys, td, psi, par, prop_psi, prop_chi, j1)
296 stop_loop = iteration_manager(sys%namespace, j1, par, par_prev, iterator)
297 if (clean_stop(sys%mc%master_comm) .or. stop_loop) exit ctr_loop
298 end do ctr_loop
299
300 call opt_control_state_end(psi)
301 call oct_prop_end(prop_chi)
302 call oct_prop_end(prop_psi)
303 call controlfunction_end(par_new)
304 call controlfunction_end(par_prev)
306 end subroutine scheme_mt03
307 ! ---------------------------------------------------------
308
309
310 ! ---------------------------------------------------------
311 subroutine scheme_wg05()
312 type(opt_control_state_t) :: psi
314
315 call oct_prop_init(prop_chi, sys%namespace, "chi", sys%gr, sys%mc)
316 call oct_prop_init(prop_psi, sys%namespace, "psi", sys%gr, sys%mc)
317
318 if (oct%mode_fixed_fluence) then
320 end if
321
322 call opt_control_state_null(psi)
323
324 call controlfunction_copy(par_new, par)
325 ctr_loop: do
326 call controlfunction_copy(par_prev, par)
327 call f_wg05(sys, td, psi, par, prop_psi, prop_chi, j1)
328 stop_loop = iteration_manager(sys%namespace, j1, par, par_prev, iterator)
329 if (clean_stop(sys%mc%master_comm) .or. stop_loop) exit ctr_loop
330 end do ctr_loop
331
332 call opt_control_state_end(psi)
333 call oct_prop_end(prop_chi)
334 call oct_prop_end(prop_psi)
335 call controlfunction_end(par_new)
336 call controlfunction_end(par_prev)
338 end subroutine scheme_wg05
339 ! ---------------------------------------------------------
340
341
342 ! ---------------------------------------------------------
343 subroutine scheme_zbr98()
344 type(opt_control_state_t) :: qcpsi
346
347 call opt_control_state_null(qcpsi)
348 call opt_control_state_copy(qcpsi, initial_st)
349 call oct_prop_init(prop_chi, sys%namespace, "chi", sys%gr, sys%mc)
350 call oct_prop_init(prop_psi, sys%namespace, "psi", sys%gr, sys%mc)
351
352 call controlfunction_copy(par_prev, par)
353 call propagate_forward(sys, td, par, oct_target, qcpsi, prop_psi)
354 j1 = target_j1(oct_target, sys%namespace, sys%gr, sys%kpoints, qcpsi)
355 stop_loop = iteration_manager(sys%namespace, j1, par, par_prev, iterator)
356 if (clean_stop(sys%mc%master_comm) .or. stop_loop) then
357 call opt_control_state_end(qcpsi)
358 call oct_prop_end(prop_chi)
359 call oct_prop_end(prop_psi)
361 return
362 end if
363
364 call controlfunction_copy(par_new, par)
365 ctr_loop: do
366 call controlfunction_copy(par_prev, par)
367 call f_zbr98(sys, td, qcpsi, prop_psi, prop_chi, par)
368 j1 = target_j1(oct_target, sys%namespace, sys%gr, sys%kpoints, qcpsi)
369 stop_loop = iteration_manager(sys%namespace, j1, par, par_prev, iterator)
370 if (clean_stop(sys%mc%master_comm) .or. stop_loop) exit ctr_loop
371 end do ctr_loop
372
373 call opt_control_state_end(qcpsi)
374 call oct_prop_end(prop_chi)
375 call oct_prop_end(prop_psi)
377 call controlfunction_end(par_prev)
379 end subroutine scheme_zbr98
380 ! ---------------------------------------------------------
381
382
383 ! ---------------------------------------------------------
384 subroutine scheme_cg()
385 integer :: dof, ierr, maxiter
386 real(real64):: step, minvalue
387 real(real64), allocatable :: theta(:)
388 real(real64), allocatable :: x(:)
389 real(real64) :: f
390 type(opt_control_state_t) :: qcpsi
392
394
395 call opt_control_state_null(qcpsi)
396 call opt_control_state_copy(qcpsi, initial_st)
397 call propagate_forward(sys, td, par, oct_target, qcpsi)
398 f = - target_j1(oct_target, sys%namespace, sys%gr, sys%kpoints, qcpsi, sys%ions) - controlfunction_j2(par)
399 call opt_control_state_end(qcpsi)
400 call iteration_manager_direct(-f, par, iterator, sys)
401 if (oct_iterator_maxiter(iterator) == 0) then
402 ! Nothing to do.
404 return
405 end if
406
407 if (oct%random_initial_guess) call controlfunction_randomize(par)
408
409 ! Set the module pointers, so that the opt_control_cg_calc and opt_control_cg_write_info routines
410 ! can use them.
411 call controlfunction_copy(par_, par)
412 sys_ => sys
413 hm_ => sys%hm
414 td_ => td
415
416 dof = controlfunction_dof(par)
417 safe_allocate(x(1:dof))
418 safe_allocate(theta(1:dof))
419 call controlfunction_get_theta(par, theta)
420 x = theta
421
422 step = oct%direct_step * m_pi
423 maxiter = oct_iterator_maxiter(iterator) - 1
424
425 select case (oct%algorithm)
426 case (option__octscheme__oct_bfgs)
427 call minimize_multidim(minmethod_bfgs2, dof, x, step, 0.1_real64, &
428 real(oct_iterator_tolerance(iterator), real64), real(oct_iterator_tolerance(iterator), real64), &
429 maxiter, opt_control_cg_calc, opt_control_cg_write_info, minvalue, ierr)
430 case (option__octscheme__oct_cg)
431 call minimize_multidim(minmethod_fr_cg, dof, x, step, 0.1_real64, &
432 real(oct_iterator_tolerance(iterator), real64), real(oct_iterator_tolerance(iterator), real64), &
433 maxiter, opt_control_cg_calc, opt_control_cg_write_info, minvalue, ierr)
434 end select
435
436 if (ierr /= 0) then
437 if (ierr <= 1024) then
438 message(1) = "Error occurred during the GSL minimization procedure:"
439 call loct_strerror(ierr, message(2))
440 call messages_fatal(2, namespace=sys%namespace)
441 else
442 message(1) = "The optimization did not meet the convergence criterion."
443 call messages_info(1, namespace=sys%namespace)
444 end if
445 end if
446
447 call controlfunction_end(par_)
448 safe_deallocate_a(x)
449 safe_deallocate_a(theta)
451 end subroutine scheme_cg
452 ! ---------------------------------------------------------
453
454
455 ! ---------------------------------------------------------
456 subroutine scheme_direct()
457 integer :: ierr, maxiter
458 real(real64):: minvalue, step
459 real(real64), allocatable :: theta(:)
460 real(real64), allocatable :: x(:)
461 real(real64) :: f
462 integer :: dim
463 type(opt_control_state_t) :: qcpsi
464
466
468 dim = controlfunction_dof(par)
469
470 call opt_control_state_null(qcpsi)
471 call opt_control_state_copy(qcpsi, initial_st)
472 call propagate_forward(sys, td, par, oct_target, qcpsi)
473 f = - target_j1(oct_target, sys%namespace, sys%gr, sys%kpoints, qcpsi, sys%ions) - controlfunction_j2(par)
474 call opt_control_state_end(qcpsi)
475 call iteration_manager_direct(-f, par, iterator, sys)
476 if (oct_iterator_maxiter(iterator) == 0) then
477 ! Nothing to do.
479 return
480 end if
481
482 safe_allocate(x(1:dim))
483 safe_allocate(theta(1:dim))
484
485 if (oct%random_initial_guess) call controlfunction_randomize(par)
486
487 ! Set the module pointers, so that the opt_control_direct_calc and opt_control_direct_message_info routines
488 ! can use them.
489 call controlfunction_copy(par_, par)
490 sys_ => sys
491 hm_ => sys%hm
492 td_ => td
493
494 ! theta may be in single precision, whereas x is always double precision.
495 call controlfunction_get_theta(par, theta)
496 x = theta
497
498 step = oct%direct_step * m_pi
499 maxiter = oct_iterator_maxiter(iterator)
500
501 call minimize_multidim_nograd(minmethod_nmsimplex, dim, x, step, &
502 real(oct_iterator_tolerance(iterator), real64) , maxiter, &
503 opt_control_direct_calc, opt_control_direct_message_info, minvalue, ierr)
504
505 if (ierr /= 0) then
506 if (ierr <= 1024) then
507 message(1) = "Error occurred during the GSL minimization procedure:"
508 call loct_strerror(ierr, message(2))
509 call messages_fatal(2, namespace=sys%namespace)
510 else
511 message(1) = "The OCT direct optimization did not meet the convergence criterion."
512 call messages_info(1, namespace=sys%namespace)
513 end if
514 end if
515
516 call controlfunction_end(par_)
517 safe_deallocate_a(x)
518 safe_deallocate_a(theta)
520 end subroutine scheme_direct
521 ! ---------------------------------------------------------
522
523
524 ! ---------------------------------------------------------
525 subroutine scheme_nlopt()
526#if defined(HAVE_NLOPT)
527 integer :: method, dim, maxiter, ierr
528 real(real64), allocatable :: x(:), xl(:), xu(:)
529 real(real64) :: step, toldr, minimum, f
530 type(opt_control_state_t) :: qcpsi
532
534
535 call opt_control_state_null(qcpsi)
536 call opt_control_state_copy(qcpsi, initial_st)
537 call propagate_forward(sys, td, par, oct_target, qcpsi)
538 f = - target_j1(oct_target, sys%namespace, sys%gr, sys%kpoints, qcpsi, sys%ions) - controlfunction_j2(par)
539 call opt_control_state_end(qcpsi)
540 call iteration_manager_direct(-f, par, iterator, sys)
541 if (oct_iterator_maxiter(iterator) == 0) then
542 ! Nothing to do.
544 return
545 end if
546
547 dim = controlfunction_dof(par)
548 safe_allocate(x(dim))
549 safe_allocate(xl(1:dim))
550 safe_allocate(xu(1:dim))
551 call controlfunction_bounds(par, xl, xu)
552
553 if (oct%random_initial_guess) call controlfunction_randomize(par)
554
555 ! Set the module pointers, so that the calc_point and write_iter_info routines
556 ! can use them.
557 call controlfunction_copy(par_, par)
558 sys_ => sys
559 hm_ => sys%hm
560 td_ => td
561
562 call controlfunction_get_theta(par, x)
563
564 maxiter = oct_iterator_maxiter(iterator)
565 step = oct%direct_step
566 select case (oct%algorithm)
567 case (option__octscheme__oct_nlopt_bobyqa)
569 case (option__octscheme__oct_nlopt_lbfgs)
570 method = minmethod_nlopt_lbfgs
571 end select
572 toldr = oct_iterator_tolerance(iterator)
573
574 call minimize_multidim_nlopt(ierr, method, dim, x, step, toldr, maxiter, opt_control_nlopt_func, minimum, &
575 xl, xu)
576 if (ierr < 1 .or. ierr > 4) then
577 message(1) = "The nlopt minimization procedure did not find convergence, or found an error"
578 write(message(2),'(a,i5)') "Error code =", ierr
579 call messages_info(2, namespace=sys%namespace)
580 end if
581
582 call controlfunction_end(par_)
583 safe_deallocate_a(xl)
584 safe_deallocate_a(xu)
585 safe_deallocate_a(x)
587#endif
588 end subroutine scheme_nlopt
589 ! ---------------------------------------------------------
590
591 end subroutine opt_control_run_legacy
592 ! ---------------------------------------------------------
593
594
595 ! ---------------------------------------------------------
596 subroutine f_zbr98(sys, td, qcpsi, prop_psi, prop_chi, par)
597 type(electrons_t), intent(inout) :: sys
598 type(td_t), intent(inout) :: td
599 type(opt_control_state_t), intent(inout) :: qcpsi
600 type(oct_prop_t), intent(inout) :: prop_psi
601 type(oct_prop_t), intent(inout) :: prop_chi
602 type(controlfunction_t), intent(inout) :: par
603
604 type(states_elec_t) :: chi
605 type(opt_control_state_t) :: qcchi
606 type(controlfunction_t) :: par_chi
607
608 push_sub(f_zbr98)
609
610 call controlfunction_copy(par_chi, par)
611
612 call target_get_state(oct_target, chi)
613 call opt_control_state_init(qcchi, chi, sys%ions)
614 call bwd_step(sys, td, oct_target, par, par_chi, qcchi, prop_chi, prop_psi)
615 call opt_control_state_copy(qcpsi, initial_st)
616 call fwd_step(sys, td, oct_target, par, par_chi, qcpsi, prop_chi, prop_psi)
617
619 call opt_control_state_end(qcchi)
620 call controlfunction_end(par_chi)
621 pop_sub(f_zbr98)
622 end subroutine f_zbr98
623
624
625 ! ---------------------------------------------------------
626 subroutine f_wg05(sys, td, qcpsi, par, prop_psi, prop_chi, j1)
627 type(electrons_t), intent(inout) :: sys
628 type(td_t), intent(inout) :: td
629 type(opt_control_state_t), intent(inout) :: qcpsi
630 type(controlfunction_t), intent(inout) :: par
631 type(oct_prop_t), intent(inout) :: prop_psi
632 type(oct_prop_t), intent(inout) :: prop_chi
633 real(real64), intent(out) :: j1
634
635 real(real64) :: new_penalty
636 type(opt_control_state_t) :: qcchi
637 type(controlfunction_t) :: parp
638
639 push_sub(f_wg05)
640
641 if (oct_iterator_current(iterator) == 0) then
642 call opt_control_state_copy(qcpsi, initial_st)
643 call propagate_forward(sys, td, par, oct_target, qcpsi, prop_psi)
644 j1 = target_j1(oct_target, sys%namespace, sys%gr, sys%kpoints, qcpsi)
645 pop_sub(f_wg05)
646 return
647 end if
648
649 call controlfunction_copy(parp, par)
650
651 call opt_control_state_null(qcchi)
652 call opt_control_state_copy(qcchi, qcpsi)
653 call target_chi(oct_target, sys%namespace, sys%gr, sys%kpoints, qcpsi, qcchi, sys%ions)
654 call bwd_step(sys, td, oct_target, par, parp, qcchi, prop_chi, prop_psi)
655
656 call controlfunction_filter(parp, filter)
657
658 ! recalc field
659 if (oct%mode_fixed_fluence) then
661 call controlfunction_set_alpha(parp, new_penalty)
663 end if
664
665 call controlfunction_copy(par, parp)
667
668 call opt_control_state_copy(qcpsi, initial_st)
669 call propagate_forward(sys, td, par, oct_target, qcpsi, prop_psi)
670
671 j1 = target_j1(oct_target, sys%namespace, sys%gr, sys%kpoints, qcpsi)
672
673 call opt_control_state_end(qcchi)
674 call controlfunction_end(parp)
675 pop_sub(f_wg05)
676 end subroutine f_wg05
677 ! ---------------------------------------------------------
678
679
680 ! ---------------------------------------------------------
681 subroutine f_striter(sys, td, par, j1)
682 type(electrons_t), intent(inout) :: sys
683 type(td_t), intent(inout) :: td
684 type(controlfunction_t), intent(inout) :: par
685 real(real64), intent(out) :: j1
686
687 type(opt_control_state_t) :: qcpsi
688 type(opt_control_state_t) :: qcchi
689 type(controlfunction_t) :: par_chi
690 type(oct_prop_t) :: prop_chi, prop_psi
691
692 push_sub(f_striter)
693
694 call oct_prop_init(prop_chi, sys%namespace, "chi", sys%gr, sys%mc)
695 call oct_prop_init(prop_psi, sys%namespace, "psi", sys%gr, sys%mc)
696
698
699 call controlfunction_copy(par_chi, par)
700
701 ! First, a forward propagation with the input field.
702 call opt_control_state_null(qcpsi)
703 call opt_control_state_copy(qcpsi, initial_st)
704 call propagate_forward(sys, td, par, oct_target, qcpsi, prop_psi)
705
706 ! Check the performance.
707 j1 = target_j1(oct_target, sys%namespace, sys%gr, sys%kpoints, qcpsi, sys%ions)
708
709 ! Set the boundary condition for the backward propagation.
710 call opt_control_state_null(qcchi)
711 call opt_control_state_copy(qcchi, qcpsi)
712 call target_chi(oct_target, sys%namespace, sys%gr, sys%kpoints, qcpsi, qcchi, sys%ions)
713
714 ! Backward propagation, while at the same time finding the output field,
715 ! which is placed at par_chi
716 call bwd_step_2(sys, td, oct_target, par, par_chi, qcchi, prop_chi, prop_psi)
717 !if (oct%mode_fixed_fluence) call controlfunction_set_fluence(par_chi)
718
719 ! Copy par_chi to par
720 call controlfunction_copy(par, par_chi)
721
722 call opt_control_state_end(qcpsi)
723 call opt_control_state_end(qcchi)
724 call controlfunction_end(par_chi)
725 call oct_prop_end(prop_chi)
726 call oct_prop_end(prop_psi)
727
728 pop_sub(f_striter)
729 end subroutine f_striter
730 ! ---------------------------------------------------------
731
732
733 ! ---------------------------------------------------------
734 subroutine f_iter(sys, td, qcpsi, par, prop_psi, prop_chi, j1)
735 type(electrons_t), intent(inout) :: sys
736 type(td_t), intent(inout) :: td
737 type(opt_control_state_t), intent(inout) :: qcpsi
738 type(controlfunction_t), intent(inout) :: par
739 type(oct_prop_t), intent(inout) :: prop_psi
740 type(oct_prop_t), intent(inout) :: prop_chi
741 real(real64), intent(out) :: j1
742
743 type(opt_control_state_t) :: qcchi
744 type(controlfunction_t) :: par_chi
745
746 push_sub(f_iter)
747
748 if (oct_iterator_current(iterator) == 0) then
749 call opt_control_state_copy(qcpsi, initial_st)
750 call propagate_forward(sys, td, par, oct_target, qcpsi, prop_psi)
751 j1 = target_j1(oct_target, sys%namespace, sys%gr, sys%kpoints, qcpsi)
752 pop_sub(f_iter)
753 return
754 end if
755
756 call controlfunction_copy(par_chi, par)
757
758 call opt_control_state_null(qcchi)
759 call opt_control_state_copy(qcchi, qcpsi)
760 call target_chi(oct_target, sys%namespace, sys%gr, sys%kpoints, qcpsi, qcchi, sys%ions)
761 call bwd_step(sys, td, oct_target, par, par_chi, qcchi, prop_chi, prop_psi)
762
763 call opt_control_state_copy(qcpsi, initial_st)
764 call fwd_step(sys, td, oct_target, par, par_chi, qcpsi, prop_chi, prop_psi)
765
766 j1 = target_j1(oct_target, sys%namespace, sys%gr, sys%kpoints, qcpsi)
767
768 call opt_control_state_end(qcchi)
769 call controlfunction_end(par_chi)
770 pop_sub(f_iter)
771 end subroutine f_iter
772 ! ---------------------------------------------------------
773
774#include "opt_control_c_inc.F90"
775#include "check_input_inc.F90"
776#include "finalcheck_inc.F90"
777
778
779end module opt_control_oct_m
780
781!! Local Variables:
782!! mode: f90
783!! coding: utf-8
784!! End:
Define which routines can be seen from the outside.
Definition: loct.F90:147
double sqrt(double __x) __attribute__((__nothrow__
subroutine minimize_multidim(method, dim, x, step, line_tol, tolgrad, toldr, maxiter, f, write_iter_info, minimum, ierr)
Definition: minimizer.F90:403
This module contains the definition of the data type that holds a "control function" used for OCT run...
subroutine, public controlfunction_filter(par, filter)
subroutine, public controlfunction_apply_envelope(cp)
subroutine, public controlfunction_set_fluence(par)
subroutine, public controlfunction_bounds(par, lower_bounds, upper_bounds)
subroutine, public controlfunction_to_realtime(par)
integer pure function, public controlfunction_dof(par)
real(real64) function, public controlfunction_fluence(par)
subroutine, public controlfunction_set_alpha(par, alpha)
subroutine, public controlfunction_write(filename, cp, namespace)
subroutine, public controlfunction_prepare_initial(par)
"Prepares" the initial guess control field: maybe it has to be normalized to a certain fluence,...
real(real64) pure function, public controlfunction_w0(par)
subroutine, public controlfunction_end(cp)
subroutine, public controlfunction_get_theta(par, theta)
real(real64) function, public controlfunction_j2(par)
subroutine, public controlfunction_randomize(par)
real(real64) pure function, public controlfunction_alpha(par, ipar)
real(real64) pure function, public controlfunction_targetfluence()
subroutine, public controlfunction_copy(cp_out, cp_in)
subroutine, public controlfunction_set_rep(par)
Transforms the control function to frequency space, if this is the space in which the functions are d...
subroutine, public controlfunction_mod_close()
subroutine, public controlfunction_to_h(cp, ext_partners)
subroutine, public controlfunction_init(cp, dt, ntiter)
Before using an controlfunction_t variable, it needs to be initialized, either by calling controlfunc...
subroutine, public controlfunction_mod_init(ext_partners, namespace, dt, max_iter, mode_fixed_fluence)
Initializes the module, should be the first subroutine to be called (the last one should be controlfu...
subroutine, public controlfunction_set(cp, ext_partners)
The external fields defined in epot_t "ep" are transferred to the control functions described in "cp"...
type(lasers_t) function, pointer, public list_get_lasers(partners)
subroutine, public filter_init(steps, namespace, dt, filter)
Definition: filter.F90:154
subroutine, public filter_write(filter, namespace)
Definition: filter.F90:333
subroutine, public filter_end(filter)
Definition: filter.F90:367
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:185
character(len= *), parameter, public oct_dir
Definition: global.F90:258
This module implements the underlying real-space grid.
Definition: grid.F90:117
integer, parameter, public independent_particles
subroutine, public initial_state_init(sys, qcstate)
Definition: initst.F90:154
Definition: io.F90:114
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:354
subroutine, public laser_write_info(lasers, namespace, dt, max_iter, iunit)
Definition: lasers.F90:926
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:930
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1125
character(len=512), private msg
Definition: messages.F90:165
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
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
integer minmethod_fr_cg
Definition: minimizer.F90:134
integer minmethod_nmsimplex
Definition: minimizer.F90:134
integer minmethod_bfgs2
Definition: minimizer.F90:134
integer minmethod_nlopt_lbfgs
Definition: minimizer.F90:134
integer minmethod_nlopt_bobyqa
Definition: minimizer.F90:134
This module implements the basic mulsisystem class, a container system for other systems.
This module contains the definition of the oct_t data type, which contains some of the basic informat...
subroutine, public oct_read_inp(oct, namespace)
Reads, from the inp file, some global information about how the QOCT run should be.
logical function, public iteration_manager(namespace, j1, par, par_prev, iterator)
subroutine, public oct_iterator_init(iterator, namespace, par)
integer pure function, public oct_iterator_maxiter(iterator)
subroutine, public oct_iterator_end(iterator, namespace)
integer pure function, public oct_iterator_current(iterator)
subroutine, public iteration_manager_direct(j, par, iterator, sys, dx)
real(real64) pure function, public oct_iterator_tolerance(iterator)
This module contains the main procedure ("opt_control_run") that is used when optimal control runs ar...
subroutine, public opt_control_cg_calc(n, x, f, getgrad, df)
subroutine check_faulty_runmodes(sys, tr)
subroutine, public opt_control_run(system)
subroutine oct_finalcheck(sys, td)
subroutine, public opt_control_cg_write_info(iter, n, val, maxdx, maxdf, x)
interface is required by its being passed as dummy routine to minimize_multidim
subroutine f_striter(sys, td, par, j1)
subroutine, public opt_control_direct_message_info(iter, n, val, maxdx, x)
subroutine f_zbr98(sys, td, qcpsi, prop_psi, prop_chi, par)
subroutine f_iter(sys, td, qcpsi, par, prop_psi, prop_chi, j1)
subroutine opt_control_run_legacy(sys)
This is the main procedure for all types of optimal control runs. It is called from the "run" procedu...
subroutine f_wg05(sys, td, qcpsi, par, prop_psi, prop_chi, j1)
subroutine, public opt_control_function_forward(x, f)
subroutine opt_control_nlopt_func(val, n, x, grad, need_gradient, f_data)
subroutine, public opt_control_direct_calc(n, x, f)
No intents here is unfortunately required because this will be passed to newuoa routines as a dummy f...
This module holds the "opt_control_state_t" datatype, which contains a quantum-classical state.
subroutine, public opt_control_state_end(ocs)
subroutine, public opt_control_state_null(ocs)
subroutine, public opt_control_state_copy(ocsout, ocsin)
subroutine, public opt_control_state_init(ocs, qstate, ions)
subroutine, public opt_control_get_qs(qstate, ocs)
this module contains the output system
Definition: output.F90:115
subroutine, public output_states(outp, namespace, space, dir, st, gr, ions, hm, iter)
Definition: output.F90:1899
subroutine, public propagation_mod_init(niter, eta, delta, number_checkpoints, zbr98, gradients)
This subroutine must be called before any QOCT propagations are done. It simply stores in the module ...
subroutine, public oct_prop_end(prop)
subroutine, public bwd_step(sys, td, tg, par, par_chi, qcchi, prop_chi, prop_psi)
subroutine, public propagate_forward(sys, td, par, tg, qcpsi, prop, write_iter)
subroutine, public fwd_step(sys, td, tg, par, par_chi, qcpsi, prop_chi, prop_psi)
subroutine, public oct_prop_init(prop, namespace, dirname, mesh, mc)
subroutine, public bwd_step_2(sys, td, tg, par, par_chi, qcchi, prop_chi, prop_psi)
subroutine, public propagator_elec_set_scf_prop(tr, threshold)
logical function, public clean_stop(comm)
returns true if a file named stop exists
Definition: restart.F90:274
This module handles spin dimensions of the states and the k-point distribution.
subroutine, public states_elec_end(st)
finalize the states_elec_t object
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, 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
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, public target_output(tg, namespace, space, gr, dir, ions, hm, outp)
Definition: target.F90:473
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
Definition: td.F90:114
subroutine, public td_end(td)
Definition: td.F90:610
subroutine, public td_init(td, namespace, space, gr, ions, st, ks, hm, ext_partners, outp)
Definition: td.F90:247
subroutine scheme_zbr98()
subroutine scheme_cg()
subroutine scheme_direct()
subroutine scheme_straight_iteration()
subroutine scheme_mt03()
subroutine scheme_wg05()
subroutine scheme_nlopt()
This is the data type used to hold a control function.
Class describing the electron system.
Definition: electrons.F90:214
Container class for lists of system_oct_m::system_t.
This is the datatype that contains the objects that are propagated: in principle this could be both t...
The states_elec_t class contains all electronic wave functions.