Octopus
dftb.F90
Go to the documentation of this file.
1!! Copyright (C) 2020 F. Bonafé, H. Appel
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 dftb_oct_m
23 use debug_oct_m
24#ifdef HAVE_DFTBPLUS
25 use dftbplus
26#endif
27 use global_oct_m
31 use io_oct_m
33 use ions_oct_m
34 use, intrinsic :: iso_fortran_env
35 use lasers_oct_m
37 use mpi_oct_m
40 use parser_oct_m
46 use system_oct_m
48 use unit_oct_m
51
52 implicit none
53
54 private
55 public :: &
56 dftb_t, &
58
59 integer, parameter :: &
60 OUTPUT_COORDINATES = 1, &
62
65 type, extends(system_t) :: dftb_t
66 integer :: n_atom
67 real(real64), allocatable :: coords(:,:), gradients(:,:)
68 real(real64), allocatable :: acc(:,:)
69 real(real64), allocatable :: tot_force(:,:)
70 real(real64), allocatable :: vel(:,:)
71 real(real64), allocatable :: prev_tot_force(:,:)
72 integer, allocatable :: species(:)
73 integer :: dynamics
74 real(real64), allocatable :: mass(:)
75 real(real64), allocatable :: atom_charges(:,:)
76 character(len=LABEL_LEN), allocatable :: labels(:)
77 real(real64), allocatable :: prev_acc(:,:,:)
78 real(real64) :: scc_tolerance
79 class(ions_t), pointer :: ions => null()
80 type(c_ptr) :: output_handle(2)
81 type(ion_dynamics_t) :: ions_dyn
82 class(lasers_t), pointer :: lasers => null()
83 logical :: laser_field
84 real(real64) :: field(3)
85 real(real64) :: energy
86#ifdef HAVE_DFTBPLUS
87 type(TDftbPlus) :: dftbp
88#endif
89 contains
90 procedure :: init_interaction => dftb_init_interaction
91 procedure :: initialize => dftb_initialize
92 procedure :: do_algorithmic_operation => dftb_do_algorithmic_operation
93 procedure :: output_start => dftb_output_start
94 procedure :: output_write => dftb_output_write
95 procedure :: output_finish => dftb_output_finish
96 procedure :: is_tolerance_reached => dftb_is_tolerance_reached
97 procedure :: copy_quantities_to_interaction => dftb_copy_quantities_to_interaction
98 procedure :: init_interaction_as_partner => dftb_init_interaction_as_partner
99 procedure :: update_interactions_start => dftb_update_interactions_start
100 procedure :: restart_write_data => dftb_restart_write_data
101 procedure :: restart_read_data => dftb_restart_read_data
102 procedure :: update_kinetic_energy => dftb_update_kinetic_energy
103 final :: dftb_finalize
104 end type dftb_t
105
106 interface dftb_t
107 procedure dftb_constructor
108 end interface dftb_t
109
111 integer, parameter :: &
112 EHRENFEST = 1, &
113 bo = 2
114
115contains
117 ! ---------------------------------------------------------
122 function dftb_constructor(namespace, grp) result(sys)
123 class(dftb_t), pointer :: sys
124 type(namespace_t), intent(in) :: namespace
125 type(mpi_grp_t), intent(in) :: grp
126
127 push_sub(dftb_constructor)
128
129#ifndef HAVE_DFTBPLUS
130 message(1) = "DFTB+ system not available. This feature requires compiling Octopus with the DFTB+ library"
131 call messages_fatal(1, namespace=namespace)
132#endif
133
134 allocate(sys)
135
136 call dftb_init(sys, namespace, grp)
137
138 pop_sub(dftb_constructor)
139 end function dftb_constructor
140
141 ! ---------------------------------------------------------
146 ! ---------------------------------------------------------
147 subroutine dftb_init(this, namespace, grp)
148 class(dftb_t), target, intent(inout) :: this
149 type(namespace_t), intent(in) :: namespace
150 type(mpi_grp_t), intent(in) :: grp
151
152 integer :: ii, jj, ispec
153 character(len=MAX_PATH_LEN) :: slako_dir
154 character(len=1), allocatable :: max_ang_mom(:)
155 character(len=LABEL_LEN) :: this_max_ang_mom, this_label
156 integer :: n_maxang_block
157 type(block_t) :: blk
158 real(real64) :: initial_temp
159
160#ifdef HAVE_DFTBPLUS
161 type(tdftbplusinput) :: input
162 type(fnode), pointer :: proot, pgeo, pham, pdftb, pmaxang, pslakos, ptype2files, panalysis
163 type(fnode), pointer :: pparseropts
164 type(fnode), pointer :: pelecdyn, pperturb, plaser
165#endif
167 push_sub(dftb_init)
169 this%namespace = namespace
171 ! Currently this system does not support any interaction
172 allocate(this%supported_interactions(0))
173 allocate(this%supported_interactions_as_partner(0))
175 ! Quantities updated by the algorithm
176 call this%quantities%add(quantity_t("position", updated_on_demand = .false.))
177 call this%quantities%add(quantity_t("velocity", updated_on_demand = .false.))
179 call messages_print_with_emphasis(msg="DFTB+ System", namespace=namespace)
181 this%ions => ions_t(namespace, grp)
182 this%n_atom = this%ions%natoms
183 safe_allocate(this%coords(1:3, 1:this%n_atom))
184 safe_allocate(this%acc(1:3, 1:this%n_atom))
185 safe_allocate(this%vel(1:3, 1:this%n_atom))
186 safe_allocate(this%tot_force(1:3, 1:this%n_atom))
187 safe_allocate(this%prev_tot_force(1:3, 1:this%n_atom))
188 safe_allocate(this%gradients(1:3, 1:this%n_atom))
189 safe_allocate(this%species(1:this%n_atom))
190 safe_allocate(this%mass(1:this%n_atom))
191 safe_allocate(this%atom_charges(1:this%n_atom, 1))
192 safe_allocate(this%labels(1:this%ions%nspecies))
193 safe_allocate(max_ang_mom(1:this%ions%nspecies))
194 max_ang_mom = "s"
196 ispec = 1
197 this%species(1) = 1
198 this%labels(1) = trim(this%ions%atom(1)%label)
199
200 this%coords = this%ions%pos
201 ! mass is read from the default pseudopotential files
202 this%mass = this%ions%mass
203 do ii = 1, this%n_atom
204 if ((ii > 1) .and. .not. (any(this%labels(1:ispec) == this%ions%atom(ii)%label))) then
205 ispec = ispec + 1
206 this%labels(ispec) = trim(this%ions%atom(ii)%label)
207 end if
208 do jj = 1, ispec
209 if (trim(this%ions%atom(ii)%label) == trim(this%labels(jj))) then
210 this%species(ii) = jj
211 end if
212 end do
213 end do
214 this%vel = m_zero
215 this%tot_force = m_zero
216
217 !%Variable MaxAngularMomentum
218 !%Type block
219 !%Section DFTBPlusInterface
220 !%Description
221 !% Specifies the highest angular momentum for each atom type. All orbitals up
222 !% to that angular momentum will be included in the calculation.
223 !% Possible values for the angular momenta are s, p, d, f.
224 !% These are examples:
225 !%
226 !% <tt>%MaxAngularMomentum
227 !% <br>&nbsp;&nbsp;'O' | 'p'
228 !% <br>&nbsp;&nbsp;'H' | 's'
229 !% <br>%</tt>
230 !%End
231 n_maxang_block = 0
232 if (parse_block(namespace, 'MaxAngularMomentum', blk) == 0) then
233 n_maxang_block = parse_block_n(blk)
234 if (n_maxang_block /= this%ions%nspecies) then
235 call messages_input_error(namespace, "MaxAngularMomentum", "Wrong number of species.")
236 end if
237
238 do ii = 1, n_maxang_block
239 call parse_block_string(blk, ii-1, 0, this_label)
240 call parse_block_string(blk, ii-1, 1, this_max_ang_mom)
241 if (.not. any(["s","p","d","f"] == trim(adjustl(this_max_ang_mom)))) then
242 call messages_input_error(namespace, "MaxAngularMomentum", "Wrong maximum angular momentum for element"//trim(this_label))
243 end if
244 do jj = 1, this%ions%nspecies
245 if (trim(adjustl(this_label)) == trim(adjustl(this%labels(jj)))) then
246 max_ang_mom(jj) = trim(adjustl(this_max_ang_mom))
247 end if
248 end do
249 end do
250 end if
251 call parse_block_end(blk)
252
253 !%Variable SlakoDir
254 !%Type string
255 !%Default "./"
256 !%Section Execution::IO
257 !%Description
258 !% Folder containing the Slako files
259 !%End
260 call parse_variable(namespace, 'SlakoDir', './', slako_dir)
261
262
263 ! Dynamics variables
264
265 call ion_dynamics_init(this%ions_dyn, namespace, this%ions, .false.)
266
267 call parse_variable(namespace, 'TDDynamics', bo, this%dynamics)
268 call messages_print_var_option('TDDynamics', this%dynamics, namespace=namespace)
269 if (this%dynamics == bo) then
270 call ion_dynamics_unfreeze(this%ions_dyn)
271 end if
272
273 !%Variable InitialIonicTemperature
274 !%Type float
275 !%Default 0.0
276 !%Section DFTBPlusInterface
277 !%Description
278 !% If this variable is present, the ions will have initial velocities
279 !% velocities to the atoms following a Boltzmann distribution with
280 !% this temperature (in Kelvin). Used only if <tt>TDDynamics = Ehrenfest</tt>
281 !% and <tt>MoveIons = yes</tt>.
282 !%End
283 call parse_variable(namespace, 'InitialIonicTemperature', m_zero, initial_temp, unit = unit_kelvin)
284
285 this%lasers => lasers_t(namespace)
286 call lasers_parse_external_fields(this%lasers)
287 if (this%lasers%no_lasers > 0) then
288 this%laser_field = .true.
289 else
290 this%laser_field = .false.
291 end if
292
293#ifdef HAVE_DFTBPLUS
294 call tdftbplus_init(this%dftbp)
295
296 call this%dftbp%getEmptyInput(input)
297 call input%getRootNode(proot)
298 call setchild(proot, "Geometry", pgeo)
299 call setchildvalue(pgeo, "Periodic", .false.)
300 call setchildvalue(pgeo, "TypeNames", this%labels(1:this%ions%nspecies))
301 call setchildvalue(pgeo, "TypesAndCoordinates", reshape(this%species, [1, size(this%species)]), this%coords)
302 call setchild(proot, "Hamiltonian", pham)
303 call setchild(pham, "Dftb", pdftb)
304 call setchildvalue(pdftb, "Scc", .true.)
305
306 !%Variable SccTolerance
307 !%Type float
308 !%Section DFTBPlusInterface
309 !%Description
310 !% Self-consistent-charges convergence tolerance. Once this
311 !% tolerance has been achieved the SCC cycle will stop.
312 !%End
313 call parse_variable(namespace, 'SccTolerance', 1e-9_real64, this%scc_tolerance)
314 call messages_print_var_value('SccTolerance', this%scc_tolerance, namespace=namespace)
315 call setchildvalue(pdftb, "SccTolerance", this%scc_tolerance)
316
317 ! sub-block inside hamiltonian for the maximum angular momenta
318 call setchild(pdftb, "MaxAngularMomentum", pmaxang)
319 ! explicitly set the maximum angular momenta for the species
320 do ii = 1, this%ions%nspecies
321 call setchildvalue(pmaxang, this%labels(ii), max_ang_mom(ii))
322 end do
323
324 ! get the SK data
325 ! You should provide the skfiles as found in the external/slakos/origin/mio-1-1/ folder. These can
326 ! be downloaded with the utils/get_opt_externals script
327 call setchild(pdftb, "SlaterKosterFiles", pslakos)
328 call setchild(pslakos, "Type2FileNames", ptype2files)
329 call setchildvalue(ptype2files, "Prefix", slako_dir)
330 call setchildvalue(ptype2files, "Separator", "-")
331 call setchildvalue(ptype2files, "Suffix", ".skf")
332
333 ! set up analysis options
334 call setchild(proot, "Analysis", panalysis)
335 call setchildvalue(panalysis, "CalculateForces", .true.)
336
337 call setchild(proot, "ParserOptions", pparseropts)
338 call setchildvalue(pparseropts, "ParserVersion", 5)
339
340 if (this%dynamics == ehrenfest) then
341 call setchild(proot, "ElectronDynamics", pelecdyn)
342 call setchildvalue(pelecdyn, "IonDynamics", this%ions_dyn%is_active())
343 if (this%ions_dyn%is_active()) then
344 call setchildvalue(pelecdyn, "InitialTemperature", initial_temp)
345 end if
346
347 ! initialize with wrong arguments for the moment, will be overriden later
348 call setchildvalue(pelecdyn, "Steps", 1)
349 call setchildvalue(pelecdyn, "TimeStep", m_one)
350 call setchild(pelecdyn, "Perturbation", pperturb)
351 if (this%laser_field) then
352 call setchild(pperturb, "Laser", plaser)
353 call setchildvalue(plaser, "PolarizationDirection", [ m_one , m_zero , m_zero ])
354 call setchildvalue(plaser, "LaserEnergy", m_one)
355 call setchildvalue(pelecdyn, "FieldStrength", m_one)
356 else
357 call setchild(pperturb, "None", plaser)
358 end if
359 end if
360
361 message(1) = 'Input tree in HSD format:'
362 call messages_info(1, namespace=namespace)
363 call dumphsd(input%hsdTree, stdout)
364
365 ! initialise the DFTB+ calculator
366 call this%dftbp%setupCalculator(input)
367 call this%dftbp%setGeometry(this%coords)
368#endif
369
370 pop_sub(dftb_init)
371 end subroutine dftb_init
372
373 ! ---------------------------------------------------------
374 subroutine dftb_init_interaction(this, interaction)
375 class(dftb_t), target, intent(inout) :: this
376 class(interaction_t), intent(inout) :: interaction
377
378 push_sub(dftb_init_interaction)
379
380 select type (interaction)
381 class default
382 message(1) = "Trying to initialize an unsupported interaction by DFTB+."
383 call messages_fatal(1, namespace=this%namespace)
384 end select
386 pop_sub(dftb_init_interaction)
387 end subroutine dftb_init_interaction
388
389 ! ---------------------------------------------------------
390 subroutine dftb_initialize(this)
391 class(dftb_t), intent(inout) :: this
392
393 push_sub(dftb_initialize)
394
395#ifdef HAVE_DFTBPLUS
396 select type (algo => this%algo)
397 class is (propagator_t)
398 select case (this%dynamics)
399 case (bo)
400 call this%dftbp%getGradients(this%gradients)
401 this%tot_force = -this%gradients
402 case (ehrenfest)
403 call this%dftbp%getEnergy(this%energy)
404 call this%dftbp%initializeTimeProp(algo%dt, this%laser_field, .false.)
405 end select
406 end select
407#endif
408
409 pop_sub(dftb_initialize)
410 end subroutine dftb_initialize
411
412 ! ---------------------------------------------------------
413 logical function dftb_do_algorithmic_operation(this, operation, updated_quantities) result(done)
414 class(dftb_t), intent(inout) :: this
415 class(algorithmic_operation_t), intent(in) :: operation
416 character(len=:), allocatable, intent(out) :: updated_quantities(:)
417
418 integer :: ii, jj, il
419 type(tdf_t) :: ff, phi
420 complex(real64) :: amp, pol(3)
421 real(real64) :: time, omega
422
424 call profiling_in(trim(this%namespace%get())//":"//trim(operation%id))
425
426 select type (algo => this%algo)
427 class is (propagator_t)
428
429 done = .true.
430 select case (this%dynamics)
431 case (bo)
432 ! Born-Oppenheimer dynamics
433 select case (operation%id)
435 ! Do nothing
436
437 case (verlet_start)
438 safe_allocate(this%prev_acc(1:3, 1:this%n_atom, 1))
439 do jj = 1, this%n_atom
440 this%acc(1:3, jj) = this%tot_force(1:3, jj) / this%mass(jj)
441 end do
442
443 case (verlet_finish)
444 safe_deallocate_a(this%prev_acc)
445
446 case (verlet_update_pos)
447 do jj = 1, this%n_atom
448 this%coords(1:3, jj) = this%coords(1:3, jj) + algo%dt * this%vel(1:3, jj) &
449 + m_half * algo%dt**2 * this%acc(1:3, jj)
450 end do
451 updated_quantities = ["position"]
452
453 case (verlet_compute_acc)
454 do ii = size(this%prev_acc, dim=3) - 1, 1, -1
455 this%prev_acc(1:3, 1:this%n_atom, ii + 1) = this%prev_acc(1:3, 1:this%n_atom, ii)
456 end do
457 this%prev_acc(1:3, 1:this%n_atom, 1) = this%acc(1:3, 1:this%n_atom)
458#ifdef HAVE_DFTBPLUS
459 call this%dftbp%setGeometry(this%coords)
460 call this%dftbp%getGradients(this%gradients)
461 this%tot_force = -this%gradients
462#endif
463 do jj = 1, this%n_atom
464 this%acc(1:3, jj) = this%tot_force(1:3, jj) / this%mass(jj)
465 end do
466
467 case (verlet_compute_vel)
468 this%vel(1:3, 1:this%n_atom) = this%vel(1:3, 1:this%n_atom) &
469 + m_half * algo%dt * (this%prev_acc(1:3, 1:this%n_atom, 1) + &
470 this%acc(1:3, 1:this%n_atom))
471 updated_quantities = ["velocity"]
472
473 case default
474 done = .false.
475 end select
476
477 case (ehrenfest)
478 ! Ehrenfest dynamics
479 select case (operation%id)
481 ! Do nothing
482 case (verlet_start)
483 !Do nothing
484 case (verlet_finish)
485 !Do nothing
486 case (verlet_update_pos)
487 this%field = m_zero
488 time = this%iteration%value()
489 do il = 1, this%lasers%no_lasers
490 ! get properties of laser
491 call laser_get_f(this%lasers%lasers(il), ff)
492 call laser_get_phi(this%lasers%lasers(il), phi)
493 omega = laser_carrier_frequency(this%lasers%lasers(il))
494 pol = laser_polarization(this%lasers%lasers(il))
495 ! calculate electric field from laser
496 amp = tdf(ff, time) * exp(m_zi * (omega*time + tdf(phi, time)))
497 this%field(1:3) = this%field(1:3) + real(amp*pol(1:3), real64)
498 end do
499#ifdef HAVE_DFTBPLUS
500 call this%dftbp%setTdElectricField(this%field)
501 call this%dftbp%doOneTdStep(this%iteration%counter(), atomnetcharges=this%atom_charges, coord=this%coords,&
502 force=this%tot_force, energy=this%energy)
503#endif
504 case (verlet_compute_acc)
505 !Do nothing
506 case (verlet_compute_vel)
507 !Do nothing
508 case default
509 done = .false.
510 end select
511
512 end select
513
514 end select
515
516 call profiling_out(trim(this%namespace%get())//":"//trim(operation%id))
519
520 ! ---------------------------------------------------------
521 logical function dftb_is_tolerance_reached(this, tol) result(converged)
522 class(dftb_t), intent(in) :: this
523 real(real64), intent(in) :: tol
524
526
527 ! this routine is never called at present, no reason to be here
528 assert(.false.)
529 converged = .false.
530
533
534 ! ---------------------------------------------------------
535 subroutine dftb_output_start(this)
536 class(dftb_t), intent(inout) :: this
537
538 push_sub(dftb_output_start)
539
540 select type (algo => this%algo)
541 class is (propagator_t)
542 ! Create output handle
543 call io_mkdir('td.general', this%namespace)
544 if (mpi_world%is_root()) then
545 call write_iter_init(this%output_handle(output_coordinates), 0, algo%dt, &
546 trim(io_workpath("td.general/coordinates", this%namespace)))
547 call write_iter_init(this%output_handle(output_forces), 0, algo%dt, &
548 trim(io_workpath("td.general/forces", this%namespace)))
549 end if
550
551 ! Output info for first iteration
552 call this%output_write()
553 end select
554
555 pop_sub(dftb_output_start)
556 end subroutine dftb_output_start
557
558 ! ---------------------------------------------------------
559 subroutine dftb_output_finish(this)
560 class(dftb_t), intent(inout) :: this
561
562 push_sub(dftb_output_finish)
563
564 select type (algo => this%algo)
565 class is (propagator_t)
566 if (mpi_world%is_root()) then
567 call write_iter_end(this%output_handle(output_coordinates))
568 call write_iter_end(this%output_handle(output_forces))
569 end if
570 end select
571
572 pop_sub(dftb_output_finish)
573 end subroutine dftb_output_finish
574
575 ! ---------------------------------------------------------
576 subroutine dftb_output_write(this)
577 class(dftb_t), intent(inout) :: this
578
579 integer :: idir, iat, iout
580 character(len=50) :: aux
581 character(1) :: out_label(2)
582 real(real64) :: tmp(3)
583
584 if (.not. mpi_world%is_root()) return ! only first node outputs
585
586 push_sub(dftb_output_write)
587 call profiling_in(trim(this%namespace%get())//":"//"OUTPUT_WRITE")
588
589 select type (algo => this%algo)
590 class is (propagator_t)
591 out_label(1) = "x"
592 out_label(2) = "f"
593
594 if (this%iteration%counter() == 0) then
595 ! header
596 do iout = 1, 2
597 call write_iter_clear(this%output_handle(iout))
598 call write_iter_string(this%output_handle(iout),'#####################################################################')
599 call write_iter_nl(this%output_handle(iout))
600 call write_iter_string(this%output_handle(iout),'# HEADER')
601 call write_iter_nl(this%output_handle(iout))
602
603 ! first line: column names
604 call write_iter_header_start(this%output_handle(iout))
605
606 do iat = 1, this%n_atom
607 do idir = 1, 3
608 write(aux, '(a1,a1,i3,a1,i3,a1)') out_label(iout),'(', iat, ',', idir, ')'
609 call write_iter_header(this%output_handle(iout), aux)
610 end do
611 end do
612 call write_iter_nl(this%output_handle(iout))
613
614 ! second line: units
615 call write_iter_string(this%output_handle(iout), '#[Iter n.]')
616 call write_iter_header(this%output_handle(iout), '[' // trim(units_abbrev(units_out%time)) // ']')
617 end do
618
619 do iat = 1, this%n_atom
620 do idir = 1, 3
621 call write_iter_header(this%output_handle(output_coordinates), '[' // trim(units_abbrev(units_out%length)) // ']')
622 call write_iter_header(this%output_handle(output_forces), '[' // trim(units_abbrev(units_out%force)) // ']')
623 end do
624 end do
625
626 do iout = 1, 2
627 call write_iter_nl(this%output_handle(iout))
628 call write_iter_string(this%output_handle(iout),'#######################################################################')
629 call write_iter_nl(this%output_handle(iout))
630 end do
631 end if
632
633 call write_iter_start(this%output_handle(output_coordinates))
634 call write_iter_start(this%output_handle(output_forces))
635
636 do iat = 1, this%n_atom
637 ! Position
638 tmp(1:3) = units_from_atomic(units_out%length, this%coords(1:3, iat))
639 call write_iter_double(this%output_handle(output_coordinates), tmp, 3)
640 ! Force
641 tmp(1:3) = units_from_atomic(units_out%force, this%tot_force(1:3, iat))
642 call write_iter_double(this%output_handle(output_forces), tmp, 3)
643 end do
644
645 call write_iter_nl(this%output_handle(output_coordinates))
646 call write_iter_nl(this%output_handle(output_forces))
647
648 end select
649
650 call profiling_out(trim(this%namespace%get())//":"//"OUTPUT_WRITE")
651 pop_sub(dftb_output_write)
652 end subroutine dftb_output_write
653
654 ! ---------------------------------------------------------
655 subroutine dftb_init_interaction_as_partner(partner, interaction)
656 class(dftb_t), intent(in) :: partner
657 class(interaction_surrogate_t), intent(inout) :: interaction
658
660
661 select type (interaction)
662 class default
663 message(1) = "Unsupported interaction."
664 call messages_fatal(1, namespace=partner%namespace)
665 end select
666
669
670 ! ---------------------------------------------------------
671 subroutine dftb_copy_quantities_to_interaction(partner, interaction)
672 class(dftb_t), intent(inout) :: partner
673 class(interaction_surrogate_t), intent(inout) :: interaction
674
676
677 select type (interaction)
678 class default
679 message(1) = "Unsupported interaction."
680 call messages_fatal(1, namespace=partner%namespace)
681 end select
682
685
686 ! ---------------------------------------------------------
687 subroutine dftb_update_interactions_start(this)
688 class(dftb_t), intent(inout) :: this
689
691
692 ! Store previous force, as it is used as SCF criterium
693 this%prev_tot_force(1:3, 1:this%n_atom) = this%tot_force(1:3, 1:this%n_atom)
694
697
698 ! ---------------------------------------------------------
699 subroutine dftb_update_kinetic_energy(this)
700 class(dftb_t), intent(inout) :: this
701
703
704 this%kinetic_energy = m_zero
705
708
709 ! ---------------------------------------------------------
710 subroutine dftb_restart_write_data(this)
711 class(dftb_t), intent(inout) :: this
712
714
715 message(1) = "DFTB system "//trim(this%namespace%get())//" cannot write restart data."
716 call messages_warning(1, namespace=this%namespace)
717
719 end subroutine dftb_restart_write_data
721 ! ---------------------------------------------------------
722 ! this function returns true if restart data could be read
723 logical function dftb_restart_read_data(this)
724 class(dftb_t), intent(inout) :: this
725
726 push_sub(dftb_restart_read_data)
727
728 ! restarting not yet supported
729 dftb_restart_read_data = .false.
730
733
734 ! ---------------------------------------------------------
735 subroutine dftb_finalize(this)
736 type(dftb_t), intent(inout) :: this
737
738 push_sub(dftb_finalize)
739
740 safe_deallocate_a(this%coords)
741 safe_deallocate_a(this%acc)
742 safe_deallocate_a(this%vel)
743 safe_deallocate_a(this%tot_force)
744 safe_deallocate_a(this%prev_tot_force)
745 safe_deallocate_a(this%gradients)
746 safe_deallocate_a(this%species)
747 safe_deallocate_a(this%mass)
748 call ion_dynamics_end(this%ions_dyn)
749
750 deallocate(this%ions)
751
752 ! No call to safe_deallocate macro here, as it gives an ICE with gfortran
753 if (associated(this%lasers)) then
754 deallocate(this%lasers)
755 end if
756
757#ifdef HAVE_DFTBPLUS
758 call tdftbplus_destruct(this%dftbp)
759#endif
760
761 call system_end(this)
762
763 pop_sub(dftb_finalize)
764 end subroutine dftb_finalize
765
766end module dftb_oct_m
767
768!! Local Variables:
769!! mode: f90
770!! coding: utf-8
771!! End:
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:182
Writes to the corresponding file and adds one to the iteration. Must be called after write_iter_init(...
Definition: write_iter.F90:163
double exp(double __x) __attribute__((__nothrow__
This module implements the basic elements defining algorithms.
Definition: algorithm.F90:143
subroutine dftb_output_write(this)
Definition: dftb.F90:574
subroutine dftb_update_interactions_start(this)
Definition: dftb.F90:685
subroutine dftb_init_interaction_as_partner(partner, interaction)
Definition: dftb.F90:653
subroutine dftb_output_finish(this)
Definition: dftb.F90:557
subroutine dftb_restart_write_data(this)
Definition: dftb.F90:708
subroutine dftb_init_interaction(this, interaction)
Definition: dftb.F90:386
logical function dftb_restart_read_data(this)
Definition: dftb.F90:721
logical function dftb_do_algorithmic_operation(this, operation, updated_quantities)
Definition: dftb.F90:411
subroutine dftb_update_kinetic_energy(this)
Definition: dftb.F90:697
subroutine dftb_initialize(this)
Definition: dftb.F90:402
integer, parameter output_forces
Definition: dftb.F90:154
class(dftb_t) function, pointer dftb_constructor(namespace, grp)
The factory routine (or constructor) allocates a pointer of the corresponding type and then calls the...
Definition: dftb.F90:218
subroutine dftb_finalize(this)
Definition: dftb.F90:733
logical function dftb_is_tolerance_reached(this, tol)
Definition: dftb.F90:519
subroutine dftb_output_start(this)
Definition: dftb.F90:533
integer, parameter bo
Definition: dftb.F90:206
subroutine, public dftb_init(this, namespace, grp)
The init routine is a module level procedure This has the advantage that different classes can have d...
Definition: dftb.F90:243
subroutine dftb_copy_quantities_to_interaction(partner, interaction)
Definition: dftb.F90:669
real(real64), parameter, public m_zero
Definition: global.F90:200
complex(real64), parameter, public m_zi
Definition: global.F90:214
real(real64), parameter, public m_half
Definition: global.F90:206
real(real64), parameter, public m_one
Definition: global.F90:201
This module defines the abstract interaction_t class, and some auxiliary classes for interactions.
Definition: io.F90:116
character(len=max_path_len) function, public io_workpath(path, namespace)
construct path name from given name and namespace
Definition: io.F90:318
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:361
subroutine, public ion_dynamics_unfreeze(this)
Unfreezes the ionic movement.
subroutine, public ion_dynamics_init(this, namespace, ions, symmetrize, symm)
subroutine, public ion_dynamics_end(this)
complex(real64) function, dimension(3), public laser_polarization(laser)
Definition: lasers.F90:730
subroutine, public lasers_parse_external_fields(this)
Definition: lasers.F90:246
subroutine, public laser_get_f(laser, ff)
Definition: lasers.F90:773
real(real64) function, public laser_carrier_frequency(laser)
Definition: lasers.F90:623
subroutine, public laser_get_phi(laser, phi)
Definition: lasers.F90:801
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:898
character(len=512), private msg
Definition: messages.F90:167
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:525
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:691
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:272
This module implements the multisystem debug functionality.
subroutine, public parse_block_string(blk, l, c, res, convert_to_c)
Definition: parser.F90:818
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:623
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
This module implements the basic propagator framework.
Definition: propagator.F90:119
character(len=algo_label_len), parameter, public store_current_status
Definition: propagator.F90:176
character(len=30), parameter, public verlet_start
character(len=30), parameter, public verlet_compute_acc
character(len=30), parameter, public verlet_update_pos
character(len=30), parameter, public verlet_finish
character(len=30), parameter, public verlet_compute_vel
This module defines the quantity_t class and the IDs for quantities, which can be exposed by a system...
Definition: quantity.F90:140
This module implements the abstract system type.
Definition: system.F90:120
subroutine, public system_end(this)
Definition: system.F90:1151
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:225
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
type(unit_t), public unit_kelvin
For converting energies into temperatures.
Explicit interfaces to C functions, defined in write_iter_low.cc.
Definition: write_iter.F90:116
subroutine, public write_iter_header(out, string)
Definition: write_iter.F90:249
subroutine, public write_iter_string(out, string)
Definition: write_iter.F90:265
subroutine, public write_iter_init(out, iter, factor, file)
Definition: write_iter.F90:229
Descriptor of one algorithmic operation.
Definition: algorithm.F90:165
class for a tight binding
Definition: dftb.F90:160
abstract interaction class
surrogate interaction class to avoid circular dependencies between modules.
Abstract class implementing propagators.
Definition: propagator.F90:144
Systems (system_t) can expose quantities that can be used to calculate interactions with other system...
Definition: quantity.F90:173
Abstract class for systems.
Definition: system.F90:174
int true(void)
subroutine input()
Definition: vdw.F90:257