69#include "../base/base_uses.f90"
88 SUBROUTINE cp_eval_at(gopt_env, x, f, gradient, master, &
89 final_evaluation, para_env)
95 TYPE(gopt_f_type),
POINTER :: gopt_env
96 REAL(KIND=
dp),
DIMENSION(:),
POINTER :: x
97 REAL(KIND=
dp),
INTENT(out),
OPTIONAL :: f
98 REAL(KIND=
dp),
DIMENSION(:),
OPTIONAL, &
100 INTEGER,
INTENT(IN) :: master
101 LOGICAL,
INTENT(IN),
OPTIONAL :: final_evaluation
102 TYPE(mp_para_env_type),
POINTER :: para_env
108 LOGICAL,
PRIVATE,
PARAMETER :: debug_this_module = .true.
109 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'gopt_f_methods'
128 TYPE(gopt_f_type),
POINTER :: gopt_env
129 REAL(kind=dp),
DIMENSION(:),
POINTER :: x0
131 INTEGER :: i, idg, j, nparticle
138 SELECT CASE (gopt_env%type_id)
142 IF (gopt_env%force_env%in_use ==
use_qmmm) &
144 IF (gopt_env%force_env%in_use ==
use_qmmmx) &
147 ALLOCATE (x0(3*nparticle))
150 SELECT CASE (gopt_env%cell_method_id)
152 CALL force_env_get(gopt_env%force_env, subsys=subsys, cell=cell)
154 gopt_env%h_ref = cell%hmat
156 IF (gopt_env%force_env%in_use ==
use_qmmm) &
158 IF (gopt_env%force_env%in_use ==
use_qmmmx) &
161 ALLOCATE (x0(3*nparticle + 6))
167 x0(idg) = cell%hmat(j, i)
177 x0(idg) = cell%hmat(j, i)
197 INTEGER,
INTENT(IN) :: its, output_unit
199 IF (output_unit > 0)
THEN
200 WRITE (unit=output_unit, fmt=
"(/,T2,26('-'))")
201 WRITE (unit=output_unit, fmt=
"(T2,A,I6)")
"OPTIMIZATION STEP: ", its
202 WRITE (unit=output_unit, fmt=
"(T2,26('-'))")
218 SUBROUTINE gopt_f_io_init(gopt_env, output_unit, opt_energy, wildcard, its, used_time)
219 TYPE(gopt_f_type),
POINTER :: gopt_env
220 INTEGER,
INTENT(IN) :: output_unit
221 REAL(kind=dp) :: opt_energy
222 CHARACTER(LEN=5) :: wildcard
223 INTEGER,
INTENT(IN) :: its
224 REAL(kind=dp) :: used_time
226 REAL(kind=dp) :: pres_int
228 SELECT CASE (gopt_env%type_id)
231 IF (.NOT. gopt_env%dimer_rotation)
THEN
232 CALL write_cycle_infos(output_unit, it=its, etot=opt_energy, wildcard=wildcard, &
235 CALL write_rot_cycle_infos(output_unit, it=its, etot=opt_energy, dimer_env=gopt_env%dimer_env, &
236 wildcard=wildcard, used_time=used_time)
240 pres_int = gopt_env%cell_env%pres_int
241 CALL write_cycle_infos(output_unit, it=its, etot=opt_energy, pres_int=pres_int, wildcard=wildcard, &
244 CALL write_cycle_infos(output_unit, it=its, etot=opt_energy, wildcard=wildcard, &
273 SUBROUTINE gopt_f_io(gopt_env, force_env, root_section, its, opt_energy, &
274 output_unit, eold, emin, wildcard, gopt_param, ndf, dx, xi, conv, pred, rat, &
275 step, rad, used_time)
276 TYPE(gopt_f_type),
POINTER :: gopt_env
279 INTEGER,
INTENT(IN) :: its
280 REAL(kind=dp),
INTENT(IN) :: opt_energy
281 INTEGER,
INTENT(IN) :: output_unit
282 REAL(kind=dp) :: eold, emin
283 CHARACTER(LEN=5) :: wildcard
285 INTEGER,
INTENT(IN),
OPTIONAL :: ndf
286 REAL(kind=dp),
DIMENSION(:),
INTENT(IN),
OPTIONAL :: dx
287 REAL(kind=dp),
DIMENSION(:),
OPTIONAL,
POINTER :: xi
288 LOGICAL,
OPTIONAL :: conv
289 REAL(kind=dp),
INTENT(IN),
OPTIONAL :: pred, rat, step, rad
290 REAL(kind=dp) :: used_time
292 INTEGER(KIND=int_8) :: max_memory
293 LOGICAL :: print_memory
294 REAL(kind=dp) :: pres_diff, pres_diff_constr, pres_int, &
296 TYPE(mp_para_env_type),
POINTER :: para_env
301 IF (print_memory)
THEN
306 SELECT CASE (gopt_env%type_id)
309 IF (.NOT. gopt_env%dimer_rotation)
THEN
310 CALL geo_opt_io(force_env=force_env, root_section=root_section, &
311 motion_section=gopt_env%motion_section, its=its, opt_energy=opt_energy)
312 CALL write_cycle_infos(output_unit, its, etot=opt_energy, ediff=opt_energy - eold, &
313 pred=pred, rat=rat, step=step, rad=rad, emin=emin, wildcard=wildcard, used_time=used_time)
315 IF (
PRESENT(conv))
THEN
316 cpassert(
PRESENT(ndf))
317 cpassert(
PRESENT(dx))
318 cpassert(
PRESENT(xi))
319 CALL check_converg(ndf, dx, xi, output_unit, conv, gopt_param, max_memory)
323 CALL write_restart(force_env=force_env, root_section=root_section)
324 CALL write_rot_cycle_infos(output_unit, its, opt_energy, opt_energy - eold, emin, gopt_env%dimer_env, &
325 used_time=used_time, wildcard=wildcard)
327 IF (
PRESENT(conv))
THEN
328 cpassert(
ASSOCIATED(gopt_env%dimer_env))
329 CALL check_rot_conv(gopt_env%dimer_env, output_unit, conv)
334 pres_diff = gopt_env%cell_env%pres_int - gopt_env%cell_env%pres_ext
335 pres_int = gopt_env%cell_env%pres_int
336 pres_tol = gopt_env%cell_env%pres_tol
337 CALL geo_opt_io(force_env=force_env, root_section=root_section, &
338 motion_section=gopt_env%motion_section, its=its, opt_energy=opt_energy)
339 CALL write_cycle_infos(output_unit, its, etot=opt_energy, ediff=opt_energy - eold, &
340 pred=pred, rat=rat, step=step, rad=rad, emin=emin, pres_int=pres_int, wildcard=wildcard, &
343 IF (
PRESENT(conv))
THEN
344 cpassert(
PRESENT(ndf))
345 cpassert(
PRESENT(dx))
346 cpassert(
PRESENT(xi))
347 IF (gopt_env%cell_env%constraint_id ==
fix_none)
THEN
348 CALL check_converg(ndf, dx, xi, output_unit, conv, gopt_param, max_memory, pres_diff, pres_tol)
350 pres_diff_constr = gopt_env%cell_env%pres_constr - gopt_env%cell_env%pres_ext
351 CALL check_converg(ndf, dx, xi, output_unit, conv, gopt_param, max_memory, pres_diff, &
352 pres_tol, pres_diff_constr)
356 CALL write_cycle_infos(output_unit, its, etot=opt_energy, ediff=opt_energy - eold, &
357 pred=pred, rat=rat, step=step, rad=rad, emin=emin, wildcard=wildcard, used_time=used_time)
359 IF (
PRESENT(conv))
THEN
360 cpassert(
PRESENT(ndf))
361 cpassert(
PRESENT(dx))
362 cpassert(
PRESENT(xi))
363 CALL check_converg(ndf, dx, xi, output_unit, conv, gopt_param, max_memory)
383 para_env, master, output_unit)
384 TYPE(gopt_f_type),
POINTER :: gopt_env
386 REAL(kind=dp),
DIMENSION(:),
POINTER :: x0
390 TYPE(mp_para_env_type),
POINTER :: para_env
391 INTEGER,
INTENT(IN) :: master, output_unit
393 IF (gopt_env%eval_opt_geo)
THEN
394 IF (.NOT. gopt_env%dimer_rotation)
THEN
395 CALL write_final_info(output_unit, conv, its, gopt_env, x0, master, &
396 para_env, force_env, gopt_env%motion_section, root_section)
399 CALL write_restart(force_env=force_env, root_section=root_section)
420 SUBROUTINE write_cycle_infos(output_unit, it, etot, ediff, pred, rat, step, rad, emin, &
421 pres_int, wildcard, used_time)
423 INTEGER,
INTENT(IN) :: output_unit, it
424 REAL(kind=dp),
INTENT(IN) :: etot
425 REAL(kind=dp),
INTENT(IN),
OPTIONAL :: ediff, pred, rat, step, rad, emin, &
427 CHARACTER(LEN=5),
OPTIONAL :: wildcard
428 REAL(kind=dp),
INTENT(IN) :: used_time
430 REAL(kind=dp) :: tmp_r1
432 IF (output_unit > 0)
THEN
433 WRITE (unit=output_unit, fmt=
"(/,T2,8('-'),A,I5,1X,12('-'))") &
434 " Informations at step = ", it
435 WRITE (unit=output_unit, fmt=
"(T2,A,T47,A)") &
436 " Optimization Method = ", wildcard
437 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
438 " Total Energy = ", etot
439 IF (
PRESENT(pres_int))
THEN
441 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
442 " Internal Pressure [bar] = ", tmp_r1
444 IF (
PRESENT(ediff))
THEN
445 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
446 " Real energy change = ", ediff
448 IF (
PRESENT(pred))
THEN
449 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
450 " Predicted change in energy = ", pred
452 IF (
PRESENT(rat))
THEN
453 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
454 " Scaling factor = ", rat
456 IF (
PRESENT(step))
THEN
457 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
458 " Step size = ", step
460 IF (
PRESENT(rad))
THEN
461 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
462 " Trust radius = ", rad
464 IF (
PRESENT(emin))
THEN
465 IF (etot < emin)
THEN
466 WRITE (unit=output_unit, fmt=
"(T2,2A)") &
467 " Decrease in energy = ", &
470 WRITE (unit=output_unit, fmt=
"(T2,2A)") &
471 " Decrease in energy = ", &
475 WRITE (unit=output_unit, fmt=
"(T2,A,F20.3)") &
476 " Used time = ", used_time
477 IF (it == 0)
WRITE (unit=output_unit, fmt=
"(T2,51('-'))")
479 END SUBROUTINE write_cycle_infos
494 SUBROUTINE write_rot_cycle_infos(output_unit, it, etot, ediff, emin, dimer_env, used_time, wildcard)
496 INTEGER,
INTENT(IN) :: output_unit, it
497 REAL(kind=dp),
INTENT(IN) :: etot
498 REAL(kind=dp),
INTENT(IN),
OPTIONAL :: ediff, emin
500 REAL(kind=dp) :: used_time
501 CHARACTER(LEN=5),
OPTIONAL :: wildcard
503 IF (output_unit > 0)
THEN
504 WRITE (unit=output_unit, fmt=
"(/,T2,4('-'),A,I5,1X,5('-'))") &
505 " Informations at rotational step = ", it
506 WRITE (unit=output_unit, fmt=
"(T2,A,T47,A)") &
507 " Optimization Method = ", wildcard
508 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
509 " Local Curvature = ", dimer_env%rot%curvature
510 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
511 " Total Rotational Force = ", etot
512 IF (
PRESENT(ediff))
THEN
513 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
514 " Real Force change = ", ediff
516 IF (
PRESENT(emin))
THEN
517 IF (etot < emin)
THEN
518 WRITE (unit=output_unit, fmt=
"(T2,2A)") &
519 " Decrease in rotational force = ", &
522 WRITE (unit=output_unit, fmt=
"(T2,2A)") &
523 " Decrease in rotational force = ", &
527 WRITE (unit=output_unit, fmt=
"(T2,A,F20.3)") &
528 " Used time = ", used_time
529 IF (it == 0)
WRITE (unit=output_unit, fmt=
"(T2,51('-'))")
531 END SUBROUTINE write_rot_cycle_infos
546 SUBROUTINE check_converg(ndf, dr, g, output_unit, conv, gopt_param, max_memory, pres_diff, &
547 pres_tol, pres_diff_constr)
549 INTEGER,
INTENT(IN) :: ndf
550 REAL(kind=dp),
INTENT(IN) :: dr(ndf), g(ndf)
551 INTEGER,
INTENT(IN) :: output_unit
552 LOGICAL,
INTENT(OUT) :: conv
554 INTEGER(KIND=int_8),
INTENT(IN) :: max_memory
555 REAL(kind=dp),
INTENT(IN),
OPTIONAL :: pres_diff, pres_tol, pres_diff_constr
558 LOGICAL :: conv_dx, conv_g, conv_p, conv_rdx, &
560 REAL(kind=dp) :: dumm, dxcon, gcon, maxdum(4), rmsgcon, &
563 dxcon = gopt_param%max_dr
564 gcon = gopt_param%max_force
565 rmsgcon = gopt_param%rms_force
566 rmsxcon = gopt_param%rms_dr
577 IF (indf == 1) maxdum(1) = abs(dr(indf))
578 dumm = dumm + dr(indf)**2
579 IF (abs(dr(indf)) > dxcon) conv_dx = .false.
580 IF (abs(dr(indf)) > maxdum(1)) maxdum(1) = abs(dr(indf))
583 IF (dumm > (rmsxcon*rmsxcon*ndf)) conv_rdx = .false.
584 maxdum(2) = sqrt(dumm/ndf)
588 IF (indf == 1) maxdum(3) = abs(g(indf))
589 dumm = dumm + g(indf)**2
590 IF (abs(g(indf)) > gcon) conv_g = .false.
591 IF (abs(g(indf)) > maxdum(3)) maxdum(3) = abs(g(indf))
594 IF (dumm > (rmsgcon*rmsgcon*ndf)) conv_rg = .false.
595 maxdum(4) = sqrt(dumm/ndf)
597 IF (
PRESENT(pres_diff_constr) .AND.
PRESENT(pres_tol))
THEN
598 conv_p = abs(pres_diff_constr) < abs(pres_tol)
599 ELSEIF (
PRESENT(pres_diff) .AND.
PRESENT(pres_tol))
THEN
600 conv_p = abs(pres_diff) < abs(pres_tol)
603 IF (output_unit > 0)
THEN
604 WRITE (unit=output_unit, fmt=
"(/,T2,A)") &
605 " Convergence check :"
606 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
607 " Max. step size = ", maxdum(1)
608 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
609 " Conv. limit for step size = ", dxcon
611 WRITE (unit=output_unit, fmt=
"(T2,2A)") &
612 " Convergence in step size = ", &
615 WRITE (unit=output_unit, fmt=
"(T2,2A)") &
616 " Convergence in step size = ", &
619 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
620 " RMS step size = ", maxdum(2)
621 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
622 " Conv. limit for RMS step = ", rmsxcon
624 WRITE (unit=output_unit, fmt=
"(T2,2A)") &
625 " Convergence in RMS step = ", &
628 WRITE (unit=output_unit, fmt=
"(T2,2A)") &
629 " Convergence in RMS step = ", &
632 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
633 " Max. gradient = ", maxdum(3)
634 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
635 " Conv. limit for gradients = ", gcon
637 WRITE (unit=output_unit, fmt=
"(T2,2A)") &
638 " Conv. in gradients = ", &
641 WRITE (unit=output_unit, fmt=
"(T2,2A)") &
642 " Conv. for gradients = ", &
645 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
646 " RMS gradient = ", maxdum(4)
647 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
648 " Conv. limit for RMS grad. = ", rmsgcon
650 WRITE (unit=output_unit, fmt=
"(T2,2A)") &
651 " Conv. in RMS gradients = ", &
655 WRITE (unit=output_unit, fmt=
"(T2,2A)") &
656 " Conv. for gradients = ", &
659 IF (
PRESENT(pres_diff) .AND.
PRESENT(pres_tol))
THEN
661 IF (
PRESENT(pres_diff_constr))
THEN
662 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
663 " Pressure Deviation [bar] = ", tmp_r1, &
664 " (without constraint)"
666 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
667 " Pressure Deviation [bar] = ", tmp_r1, &
670 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
671 " Pressure Deviation [bar] = ", tmp_r1
674 WRITE (unit=output_unit, fmt=
"(T2,A,F20.10)") &
675 " Pressure Tolerance [bar] = ", tmp_r1
677 WRITE (unit=output_unit, fmt=
"(T2,2A)") &
678 " Conv. for PRESSURE = ", &
681 WRITE (unit=output_unit, fmt=
"(T2,2A)") &
682 " Conv. for PRESSURE = ", &
686 WRITE (unit=output_unit, fmt=
"(T2,51('-'))")
687 IF (max_memory .NE. 0)
THEN
688 WRITE (output_unit,
'(T2,A,T73,I8)') &
689 'Estimated peak process memory after this step [MiB]', &
690 (max_memory + (1024*1024) - 1)/(1024*1024)
694 IF (conv_dx .AND. conv_rdx .AND. conv_g .AND. conv_rg .AND. conv_p) conv = .true.
696 IF ((conv) .AND. (output_unit > 0))
THEN
697 WRITE (unit=output_unit, fmt=
"(/,T2,A)") repeat(
"*", 79)
698 WRITE (unit=output_unit, fmt=
"(T2,A,T25,A,T78,A)") &
699 "***",
"GEOMETRY OPTIMIZATION COMPLETED",
"***"
700 WRITE (unit=output_unit, fmt=
"(T2,A)") repeat(
"*", 79)
703 END SUBROUTINE check_converg
713 SUBROUTINE check_rot_conv(dimer_env, output_unit, conv)
716 INTEGER,
INTENT(IN) :: output_unit
717 LOGICAL,
INTENT(OUT) :: conv
719 conv = (abs(dimer_env%rot%angle2) < dimer_env%rot%angle_tol)
720 IF (output_unit > 0)
THEN
721 WRITE (unit=output_unit, fmt=
"(/,T2,A)") &
722 " Convergence check :"
723 WRITE (unit=output_unit, fmt=
"(T2,A,F16.10)") &
724 " Predicted angle step size = ", dimer_env%rot%angle1
725 WRITE (unit=output_unit, fmt=
"(T2,A,F16.10)") &
726 " Effective angle step size = ", dimer_env%rot%angle2
727 WRITE (unit=output_unit, fmt=
"(T2,A,F16.10)") &
728 " Conv. limit for angle step size =", dimer_env%rot%angle_tol
730 WRITE (unit=output_unit, fmt=
"(T2,2A)") &
731 " Convergence in angle step size =", &
734 WRITE (unit=output_unit, fmt=
"(T2,2A)") &
735 " Convergence in angle step size =", &
738 WRITE (unit=output_unit, fmt=
"(T2,51('-'))")
740 IF ((conv) .AND. (output_unit > 0))
THEN
741 WRITE (unit=output_unit, fmt=
"(/,T2,A)") repeat(
"*", 79)
742 WRITE (unit=output_unit, fmt=
"(T2,A,T25,A,T78,A)") &
743 "***",
"ROTATION OPTIMIZATION COMPLETED",
"***"
744 WRITE (unit=output_unit, fmt=
"(T2,A)") repeat(
"*", 79)
747 END SUBROUTINE check_rot_conv
764 RECURSIVE SUBROUTINE write_final_info(output_unit, conv, it, gopt_env, x0, master, para_env, force_env, &
765 motion_section, root_section)
766 INTEGER,
INTENT(IN) :: output_unit
767 LOGICAL,
INTENT(IN) :: conv
768 INTEGER,
INTENT(INOUT) :: it
769 TYPE(gopt_f_type),
POINTER :: gopt_env
770 REAL(kind=dp),
DIMENSION(:),
POINTER :: x0
771 INTEGER,
INTENT(IN) :: master
772 TYPE(mp_para_env_type),
POINTER :: para_env
776 REAL(kind=dp) :: etot
784 particle_set => particles%els
788 CALL write_restart(force_env=force_env, root_section=root_section)
790 IF (output_unit > 0) &
791 WRITE (unit=output_unit, fmt=
"(/,T20,' Reevaluating energy at the minimum')")
793 CALL cp_eval_at(gopt_env, x0, f=etot, master=master, final_evaluation=.true., &
795 CALL write_geo_traj(force_env, root_section, it, etot)
798 END SUBROUTINE write_final_info
811 SUBROUTINE write_geo_traj(force_env, root_section, it, etot)
815 INTEGER,
INTENT(IN) :: it
816 REAL(kind=dp),
INTENT(IN) :: etot
818 LOGICAL :: shell_adiabatic, shell_present
824 NULLIFY (atomic_kinds)
825 NULLIFY (atomic_kind_set)
826 NULLIFY (core_particles)
827 NULLIFY (shell_particles)
832 CALL write_trajectory(force_env, root_section, it, 0.0_dp, 0.0_dp, etot,
"FORCES", middle_name=
"frc")
835 atomic_kind_set => atomic_kinds%els
837 shell_present=shell_present, &
838 shell_adiabatic=shell_adiabatic)
839 IF (shell_present)
THEN
841 core_particles=core_particles, &
842 shell_particles=shell_particles)
843 CALL write_trajectory(force_env, root_section, it=it, time=0.0_dp, dtime=0.0_dp, &
844 etot=etot, pk_name=
"SHELL_TRAJECTORY", middle_name=
"shpos", &
845 particles=shell_particles)
846 IF (shell_adiabatic)
THEN
847 CALL write_trajectory(force_env, root_section, it=it, time=0.0_dp, dtime=0.0_dp, &
848 etot=etot, pk_name=
"SHELL_FORCES", middle_name=
"shfrc", &
849 particles=shell_particles)
850 CALL write_trajectory(force_env, root_section, it=it, time=0.0_dp, dtime=0.0_dp, &
851 etot=etot, pk_name=
"CORE_TRAJECTORY", middle_name=
"copos", &
852 particles=core_particles)
853 CALL write_trajectory(force_env, root_section, it=it, time=0.0_dp, dtime=0.0_dp, &
854 etot=etot, pk_name=
"CORE_FORCES", middle_name=
"cofrc", &
855 particles=core_particles)
859 END SUBROUTINE write_geo_traj
871 TYPE(gopt_f_type),
POINTER :: gopt_env
872 INTEGER,
INTENT(IN) :: output_unit
873 CHARACTER(LEN=*),
INTENT(IN) :: label
875 CHARACTER(LEN=default_string_length) :: my_format, my_label
878 IF (output_unit > 0)
THEN
879 WRITE (unit=output_unit, fmt=
"(/,T2,A)") repeat(
"*", 79)
880 IF (gopt_env%dimer_rotation)
THEN
881 my_label =
"OPTIMIZING DIMER ROTATION"
883 my_label =
"STARTING "//gopt_env%tag(1:8)//
" OPTIMIZATION"
886 ix = (80 - 7 - len_trim(my_label))/2
889 WRITE (unit=output_unit, fmt=trim(my_format))
"***", trim(my_label),
"***"
891 ix = (80 - 7 - len_trim(label))/2
894 WRITE (unit=output_unit, fmt=trim(my_format))
"***", trim(label),
"***"
896 WRITE (unit=output_unit, fmt=
"(T2,A)") repeat(
"*", 79)
910 TYPE(gopt_f_type),
POINTER :: gopt_env
911 INTEGER,
INTENT(IN) :: output_unit
913 IF (output_unit > 0)
THEN
914 WRITE (unit=output_unit, fmt=
"(/,T2,A)") &
915 "*** MAXIMUM NUMBER OF OPTIMIZATION STEPS REACHED ***"
916 IF (.NOT. gopt_env%dimer_rotation)
THEN
917 WRITE (unit=output_unit, fmt=
"(T2,A)") &
918 "*** EXITING GEOMETRY OPTIMIZATION ***"
920 WRITE (unit=output_unit, fmt=
"(T2,A)") &
921 "*** EXITING ROTATION OPTIMIZATION ***"
939 SUBROUTINE geo_opt_io(force_env, root_section, motion_section, its, opt_energy)
943 INTEGER,
INTENT(IN) :: its
944 REAL(kind=dp),
INTENT(IN) :: opt_energy
951 TYPE(mp_para_env_type),
POINTER :: para_env
956 NULLIFY (para_env, atomic_kind_set, subsys, particle_set, &
957 local_particles, atomic_kinds, particles)
960 CALL write_restart(force_env=force_env, root_section=root_section)
963 CALL write_geo_traj(force_env, root_section, its, opt_energy)
966 CALL force_env_get(force_env, cell=cell, para_env=para_env, &
968 CALL cp_subsys_get(subsys=subsys, atomic_kinds=atomic_kinds, local_particles=local_particles, &
969 particles=particles, virial=virial)
970 atomic_kind_set => atomic_kinds%els
971 particle_set => particles%els
979 END SUBROUTINE geo_opt_io
993 TYPE(gopt_f_type),
POINTER :: gopt_env
995 REAL(kind=dp),
DIMENSION(:),
POINTER :: x
996 LOGICAL,
INTENT(IN) :: update_forces
998 INTEGER :: i, iatom, idg, j, natom, nparticle, &
1000 REAL(kind=dp) :: fc, fs, mass
1001 REAL(kind=dp),
DIMENSION(3) :: s
1008 NULLIFY (core_particles)
1010 NULLIFY (shell_particles)
1018 core_particles=core_particles, &
1019 particles=particles, &
1020 shell_particles=shell_particles)
1024 CALL cell_copy(cell, cell_ref, tag=
"CELL_OPT_REF")
1027 SELECT CASE (gopt_env%cell_method_id)
1030 CALL init_cell(cell_ref, hmat=gopt_env%h_ref)
1034 cpassert((
SIZE(x) == idg + 6))
1036 IF (update_forces)
THEN
1052 cell%hmat(j, i) = x(idg)
1059 SELECT CASE (gopt_env%cell_method_id)
1064 shell_index = particles%els(iatom)%shell_index
1065 IF (shell_index == 0)
THEN
1069 i = 3*(natom + shell_index - 1) + 1
1073 mass = particles%els(iatom)%atomic_kind%mass
1074 fc = core_particles%els(shell_index)%atomic_kind%shell%mass_core/mass
1075 fs = shell_particles%els(shell_index)%atomic_kind%shell%mass_shell/mass
1076 particles%els(iatom)%r(1:3) = fc*core_particles%els(shell_index)%r(1:3) + &
1077 fs*shell_particles%els(shell_index)%r(1:3)
1083 shell_index = particles%els(iatom)%shell_index
1084 IF (shell_index == 0)
THEN
1088 CALL real_to_scaled(s, core_particles%els(shell_index)%r, cell_ref)
1090 i = 3*(natom + shell_index - 1) + 1
1091 CALL real_to_scaled(s, shell_particles%els(shell_index)%r, cell_ref)
1094 mass = particles%els(iatom)%atomic_kind%mass
1095 fc = core_particles%els(shell_index)%atomic_kind%shell%mass_core/mass
1096 fs = shell_particles%els(shell_index)%atomic_kind%shell%mass_shell/mass
1097 particles%els(iatom)%r(1:3) = fc*core_particles%els(shell_index)%r(1:3) + &
1098 fs*shell_particles%els(shell_index)%r(1:3)
subroutine cp_eval_at(gopt_env, x, f, gradient, master, final_evaluation, para_env)
evaluete the potential energy and its gradients using an array with same dimension as the particle_se...
represent a simple array based list of the given type
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind_set(atomic_kind_set, atom_of_kind, kind_of, natom_of_kind, maxatom, natom, nshell, fist_potential_present, shell_present, shell_adiabatic, shell_check_distance, damping_present)
Get attributes of an atomic kind set.
Handles all functions related to the CELL.
subroutine, public init_cell(cell, hmat, periodic)
Initialise/readjust a simulation cell after hmat has been changed.
subroutine, public cell_create(cell, hmat, periodic, tag)
allocates and initializes a cell
Handles all functions related to the CELL.
subroutine, public scaled_to_real(r, s, cell)
Transform scaled cell coordinates real coordinates. r=h*s.
subroutine, public real_to_scaled(s, r, cell)
Transform real to scaled cell coordinates. s=h_inv*r.
subroutine, public cell_release(cell)
releases the given cell (see doc/ReferenceCounting.html)
subroutine, public cell_copy(cell_in, cell_out, tag)
Copy cell variable.
various routines to log and control the output. The idea is that decisions about where to log should ...
types that represent a subsys, i.e. a part of the system
subroutine, public cp_subsys_set(subsys, atomic_kinds, particles, local_particles, molecules, molecule_kinds, local_molecules, para_env, colvar_p, shell_particles, core_particles, gci, multipoles, results, cell)
sets various propreties of the subsys
subroutine, public cp_subsys_get(subsys, ref_count, atomic_kinds, atomic_kind_set, particles, particle_set, local_particles, molecules, molecule_set, molecule_kinds, molecule_kind_set, local_molecules, para_env, colvar_p, shell_particles, core_particles, gci, multipoles, natom, nparticle, ncore, nshell, nkind, atprop, virial, results, cell)
returns information about various attributes of the given subsys
subroutine, public pack_subsys_particles(subsys, f, r, s, v, fscale, cell)
Pack components of a subsystem particle sets into a single vector.
real(kind=dp) function, public cp_unit_from_cp2k(value, unit_str, defaults, power)
converts from the internal cp2k units to the given unit
Contains types used for a Dimer Method calculations.
Contains utilities for a Dimer Method calculations.
subroutine, public update_dimer_vec(dimer_env, motion_section)
Updates the orientation of the dimer vector in the input file.
stores a lists of integer that are local to a processor. The idea is that these integers represent ob...
Interface for the force calculations.
integer function, public force_env_get_natom(force_env)
returns the number of atoms
integer, parameter, public use_qmmm
integer, parameter, public use_qmmmx
recursive subroutine, public force_env_get(force_env, in_use, fist_env, qs_env, meta_env, fp_env, subsys, para_env, potential_energy, additional_potential, kinetic_energy, harmonic_shell, kinetic_shell, cell, sub_force_env, qmmm_env, qmmmx_env, eip_env, pwdft_env, globenv, input, force_env_section, method_name_id, root_section, mixed_env, nnp_env, embed_env)
returns various attributes about the force environment
integer function, public force_env_get_nparticle(force_env)
returns the number of particles in a force environment
contains a functional that calculates the energy and its derivatives for the geometry optimizer
subroutine, public print_geo_opt_header(gopt_env, output_unit, label)
...
subroutine, public gopt_f_io_init(gopt_env, output_unit, opt_energy, wildcard, its, used_time)
Handles the Output during an optimization run.
subroutine, public gopt_f_create_x0(gopt_env, x0)
returns the value of the parameters for the actual configuration
recursive subroutine, public gopt_f_io_finalize(gopt_env, force_env, x0, conv, its, root_section, para_env, master, output_unit)
Handles the Output at the end of an optimization run.
subroutine, public apply_cell_change(gopt_env, cell, x, update_forces)
Apply coordinate transformations after cell (shape) change.
subroutine, public gopt_f_io(gopt_env, force_env, root_section, its, opt_energy, output_unit, eold, emin, wildcard, gopt_param, ndf, dx, xi, conv, pred, rat, step, rad, used_time)
Handles the Output during an optimization run.
subroutine, public print_geo_opt_nc(gopt_env, output_unit)
...
subroutine, public gopt_f_ii(its, output_unit)
Prints iteration step of the optimization procedure on screen.
contains a functional that calculates the energy and its derivatives for the geometry optimizer
contains typo and related routines to handle parameters controlling the GEO_OPT module
Defines the basic variable types.
integer, parameter, public int_8
integer, parameter, public dp
integer, parameter, public default_string_length
Machine interface based on Fortran 2003 and POSIX.
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
prints all energy info per timestep to the screen or to user defined output files
integer(kind=int_8) function, public sample_memory(para_env)
Samples memory usage.
Interface to the message passing library MPI.
Output Utilities for MOTION_SECTION.
subroutine, public write_stress_tensor(virial, cell, motion_section, itimes, time, pos, act)
Prints the Stress Tensor.
subroutine, public write_simulation_cell(cell, motion_section, itimes, time, pos, act)
Prints the Simulation Cell.
subroutine, public write_trajectory(force_env, root_section, it, time, dtime, etot, pk_name, pos, act, middle_name, particles, extended_xmol_title)
Prints the information controlled by the TRAJECTORY section.
represent a simple array based list of the given type
Define methods related to particle_type.
subroutine, public write_structure_data(particle_set, cell, input_section)
Write structure data requested by a separate structure data input section to the output unit....
Define the data structure for the particle information.
subroutine, public apply_qmmm_translate(qmmm_env)
Apply translation to the full system in order to center the QM system into the QM box.
Routines used for force-mixing QM/MM calculations.
subroutine, public apply_qmmmx_translate(qmmmx_env)
Apply translation to the full system in order to center the QM system into the QM box.
subroutine, public virial_evaluate(atomic_kind_set, particle_set, local_particles, virial, igroup)
Computes the kinetic part of the pressure tensor and updates the full VIRIAL (PV)
represent a list of objects
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
represents a system: atoms, molecules, their pos,vel,...
Defines the environment for a Dimer Method calculation.
structure to store local (to a processor) ordered lists of integers.
wrapper to abstract the force evaluation of the various methods
calculates the potential energy of a system, and its derivatives
stores all the informations relevant to an mpi environment
represent a list of objects