Octopus
classical_particles.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2019 N. Tancogne-Dejean
3!! Copyright (C) 2020 M. Oliveira
4!! Copyright (C) 2021 S. Ohlmann, I. Albar, A. Obzhirov, A. Geondzhian, M. Lawan
5!! Copyright (C) 2022 M. Oliveira
6!!
7!! This program is free software; you can redistribute it and/or modify
8!! it under the terms of the GNU General Public License as published by
9!! the Free Software Foundation; either version 2, or (at your option)
10!! any later version.
11!!
12!! This program is distributed in the hope that it will be useful,
13!! but WITHOUT ANY WARRANTY; without even the implied warranty of
14!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15!! GNU General Public License for more details.
16!!
17!! You should have received a copy of the GNU General Public License
18!! along with this program; if not, write to the Free Software
19!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
20!! 02110-1301, USA.
21!!
22
23#include "global.h"
24
27 use debug_oct_m
29 use global_oct_m
32 use io_oct_m
33 use, intrinsic :: iso_fortran_env
36 use mpi_oct_m
39 use parser_oct_m
48 use space_oct_m
49 use system_oct_m
50 use unit_oct_m
52 use utils_oct_m
54
55 implicit none
56
57 private
58 public :: &
68
69 integer, parameter :: &
70 OUTPUT_COORDINATES = 1, &
72
73 type, extends(system_t), abstract :: classical_particles_t
74 private
75 type(c_ptr), public :: output_handle(2)
76 type(space_t), public :: space
77 integer, public :: np
78 real(real64), allocatable, public :: mass(:)
79 real(real64), allocatable, public :: pos(:,:)
80 real(real64), allocatable, public :: vel(:,:)
81 real(real64), allocatable, public :: tot_force(:,:)
82 real(real64), allocatable, public :: lj_epsilon(:)
83 real(real64), allocatable, public :: lj_sigma(:)
84 logical, allocatable, public :: fixed(:)
85 ! !! The default is to let the particles move.
86 type(propagator_data_t),public :: prop_data
87 contains
88 procedure :: do_algorithmic_operation => classical_particles_do_algorithmic_operation
89 procedure :: is_tolerance_reached => classical_particles_is_tolerance_reached
90 procedure :: copy_quantities_to_interaction => classical_particles_copy_quantities_to_interaction
91 procedure :: update_interactions_start => classical_particles_update_interactions_start
92 procedure :: update_interactions_finish => classical_particles_update_interactions_finish
93 procedure :: output_start => classical_particles_output_start
94 procedure :: output_write => classical_particles_output_write
95 procedure :: output_finish => classical_particles_output_finish
96 procedure :: restart_write_data => classical_particles_restart_write_data
97 procedure :: restart_read_data => classical_particles_restart_read_data
98 procedure :: update_kinetic_energy => classical_particles_update_kinetic_energy
99 procedure :: center_of_mass => classical_particles_center_of_mass
100 procedure :: center_of_mass_vel => classical_particles_center_of_mass_vel
101 procedure :: center => classical_particles_center
102 procedure :: axis_large => classical_particles_axis_large
103 procedure :: axis_inertia => classical_particles_axis_inertia
104 end type classical_particles_t
105
106
107contains
108
109 ! ---------------------------------------------------------
114 ! ---------------------------------------------------------
115 subroutine classical_particles_init(this, np)
116 class(classical_particles_t), intent(inout) :: this
117 integer, intent(in) :: np
118
121 this%np = np
122 safe_allocate(this%mass(1:np))
123 this%mass = m_zero
124 safe_allocate(this%pos(1:this%space%dim, 1:np))
125 this%pos = m_zero
126 safe_allocate(this%vel(1:this%space%dim, 1:np))
127 this%vel = m_zero
128 safe_allocate(this%tot_force(1:this%space%dim, 1:np))
129 this%tot_force = m_zero
130 safe_allocate(this%fixed(1:np))
131 this%fixed = .false. ! By default we let the particles move.
132 safe_allocate(this%lj_epsilon(1:np))
133 this%lj_epsilon = m_zero
134 safe_allocate(this%lj_sigma(1:np))
135 this%lj_sigma = m_zero
136
137 ! No interaction directly supported by the classical particles (but classes
138 ! that extend this one can add their own)
139 allocate(this%supported_interactions(0))
140 allocate(this%supported_interactions_as_partner(0))
141
142 ! Quantities updated by the algorithm
143 call this%quantities%add(quantity_t("position", updated_on_demand = .false.))
144 call this%quantities%add(quantity_t("velocity", updated_on_demand = .false.))
145
146 ! Other quantities
147 call this%quantities%add(quantity_t("mass", updated_on_demand = .false., always_available=.true.))
148
150 end subroutine classical_particles_init
151
152 ! ---------------------------------------------------------
153 subroutine classical_particles_copy(this, cp_in)
154 class(classical_particles_t), intent(out) :: this
155 class(classical_particles_t), intent(in) :: cp_in
156
158
159 this%np = cp_in%np
160 safe_allocate_source_a(this%mass, cp_in%mass)
161 safe_allocate_source_a(this%pos, cp_in%pos)
162 safe_allocate_source_a(this%vel, cp_in%vel)
163 safe_allocate_source_a(this%tot_force, cp_in%tot_force)
164 safe_allocate_source_a(this%fixed, cp_in%fixed)
165
166 this%kinetic_energy = cp_in%kinetic_energy
167
168 this%quantities = cp_in%quantities
169 this%supported_interactions = cp_in%supported_interactions
171 this%prop_data = cp_in%prop_data
176 ! ---------------------------------------------------------
177 subroutine classical_particles_init_interaction(this, interaction)
178 class(classical_particles_t), target, intent(inout) :: this
179 class(interaction_surrogate_t), intent(inout) :: interaction
180
182
183 select type (interaction)
184 class default
185 message(1) = "Trying to initialize an unsupported interaction by classical particles."
186 call messages_fatal(1, namespace=this%namespace)
187 end select
192 ! ---------------------------------------------------------
193 logical function classical_particles_do_algorithmic_operation(this, operation, updated_quantities) result(done)
194 class(classical_particles_t), intent(inout) :: this
195 class(algorithmic_operation_t), intent(in) :: operation
196 character(len=:), allocatable, intent(out) :: updated_quantities(:)
198 integer :: ii, ip
199 real(real64), allocatable :: tmp_pos(:,:,:), tmp_vel(:,:,:)
200 real(real64) :: factor
201
203 call profiling_in(trim(this%namespace%get())//":"//trim(operation%id))
204
205 done = .true.
206 select type (prop => this%algo)
207 class is (propagator_t)
208
209 select case (operation%id)
211 this%prop_data%save_pos(:, 1:this%np) = this%pos(:, 1:this%np)
212 this%prop_data%save_vel(:, 1:this%np) = this%vel(:, 1:this%np)
213
214 case (verlet_start)
215 if (.not. this%prop_data%initialized) then
216 call this%prop_data%initialize(prop, this%space%dim, this%np)
217 do ip = 1, this%np
218 if (this%fixed(ip)) then
219 this%prop_data%acc(:, ip) = m_zero
220 else
221 this%prop_data%acc(:, ip) = this%tot_force(:, ip) / this%mass(ip)
222 end if
223 end do
224 end if
225
226 case (verlet_finish)
227 call this%prop_data%end()
228
229 case (beeman_finish)
230 call this%prop_data%end()
231
232 case (verlet_update_pos)
233 this%pos(:, 1:this%np) = this%pos(:, 1:this%np) + prop%dt * this%vel(:, 1:this%np) &
234 + m_half * prop%dt**2 * this%prop_data%acc(:, 1:this%np)
235 updated_quantities = ["position"]
236
238 do ii = size(this%prop_data%prev_acc, dim=3) - 1, 1, -1
239 this%prop_data%prev_acc(:, 1:this%np, ii + 1) = this%prop_data%prev_acc(:, 1:this%np, ii)
240 end do
241 do ip = 1, this%np
242 this%prop_data%prev_acc(:, ip, 1) = this%prop_data%acc(:, ip)
243 if (this%fixed(ip)) then
244 this%prop_data%acc(:, ip) = m_zero
245 else
246 this%prop_data%acc(:, ip) = this%tot_force(:, ip) / this%mass(ip)
247 end if
248 end do
249
250 case (verlet_compute_vel)
251 this%vel(:, 1:this%np) = this%vel(:, 1:this%np) &
252 + m_half * prop%dt * (this%prop_data%prev_acc(:, 1:this%np, 1) + this%prop_data%acc(:, 1:this%np))
253 updated_quantities = ["velocity"]
254
255 case (beeman_start)
256 if (.not. this%prop_data%initialized) then
257 call this%prop_data%initialize(prop, this%space%dim, this%np)
258 do ip = 1, this%np
259 if (this%fixed(ip)) then
260 this%prop_data%acc(:, ip) = m_zero
261 else
262 this%prop_data%acc(:, ip) = this%tot_force(:, ip) / this%mass(ip)
263 end if
264 this%prop_data%prev_acc(:, ip, 1) = this%prop_data%acc(:, ip)
265 end do
266 end if
267
268 case (beeman_predict_pos)
269 this%pos(:, 1:this%np) = this%pos(:, 1:this%np) + prop%dt * this%vel(:, 1:this%np) + &
270 m_one/6.0_real64 * prop%dt**2 * (m_four*this%prop_data%acc(:, 1:this%np) - this%prop_data%prev_acc(:, 1:this%np, 1))
271
272 if (.not. prop%predictor_corrector) then
273 updated_quantities = ["position"]
274 end if
275
276 case (beeman_predict_vel)
277 this%vel(:, 1:this%np) = this%vel(:, 1:this%np) + &
278 m_one/6.0_real64 * prop%dt * ( m_two * this%prop_data%acc(:, 1:this%np) + &
279 5.0_real64 * this%prop_data%prev_acc(:, 1:this%np, 1) - this%prop_data%prev_acc(:, 1:this%np, 2))
280 updated_quantities = ["velocity"]
281
282 case (beeman_correct_pos)
283 this%pos(:, 1:this%np) = this%prop_data%save_pos(:, 1:this%np) &
284 + prop%dt * this%prop_data%save_vel(:, 1:this%np) &
285 + m_one/6.0_real64 * prop%dt**2 * (this%prop_data%acc(:, 1:this%np) &
286 + m_two * this%prop_data%prev_acc(:, 1:this%np, 1))
287 updated_quantities = ["position"]
289 case (beeman_correct_vel)
290 this%vel(:, 1:this%np) = this%prop_data%save_vel(:, 1:this%np) &
291 + m_half * prop%dt * (this%prop_data%acc(:, 1:this%np) + this%prop_data%prev_acc(:, 1:this%np, 1))
292 updated_quantities = ["velocity"]
293
294 case (expmid_2step_start)
295 if (.not. this%prop_data%initialized) then
296 call this%prop_data%initialize(prop, this%space%dim, this%np)
297 this%prop_data%prev_pos(:, 1:this%np, 1) = this%pos(:, 1:this%np)
298 this%prop_data%prev_vel(:, 1:this%np, 1) = this%vel(:, 1:this%np)
299 end if
300
302 call this%prop_data%end()
303
305 this%pos(:, 1:this%np) = 1.5_real64*this%prop_data%save_pos(:, 1:this%np) &
306 - m_half*this%prop_data%prev_pos(:, 1:this%np, 1)
307 this%vel(:, 1:this%np) = 1.5_real64*this%prop_data%save_vel(:, 1:this%np) &
308 - m_half*this%prop_data%prev_vel(:, 1:this%np, 1)
309 this%prop_data%prev_pos(:, 1:this%np, 1) = this%prop_data%save_pos(:, 1:this%np)
310 this%prop_data%prev_vel(:, 1:this%np, 1) = this%prop_data%save_vel(:, 1:this%np)
311 updated_quantities = ["position", "velocity"]
312
313 case (update_hamiltonian)
314 do ip = 1, this%np
315 if (this%fixed(ip)) then
316 this%prop_data%hamiltonian_elements(:, ip) = m_zero
317 else
318 this%prop_data%hamiltonian_elements(:, ip) = this%tot_force(:, ip) / (this%mass(ip) * this%pos(:, ip) + m_epsilon)
319 end if
320 end do
321
323 safe_allocate(tmp_pos(1:this%space%dim, 1:this%np, 1:2))
324 safe_allocate(tmp_vel(1:this%space%dim, 1:this%np, 1:2))
325 ! apply exponential - at some point this could use the machinery of
326 ! exponential_apply (but this would require a lot of boilerplate code
327 ! like a Hamiltonian class etc)
328 ! prop_data%save_pos/vel contain the state at t - this is the state we want to
329 ! apply the Hamiltonian to
330 tmp_pos(:, 1:this%np, 1) = this%prop_data%save_pos(:, 1:this%np)
331 tmp_vel(:, 1:this%np, 1) = this%prop_data%save_vel(:, 1:this%np)
332 this%pos(:, 1:this%np) = this%prop_data%save_pos(:, 1:this%np)
333 this%vel(:, 1:this%np) = this%prop_data%save_vel(:, 1:this%np)
334 ! compute exponential with Taylor expansion
335 factor = m_one
336 do ii = 1, 4
337 factor = factor * prop%dt / ii
338 do ip = 1, this%np
339 ! apply hamiltonian
340 tmp_pos(:, ip, 2) = tmp_vel(:, ip, 1)
341 tmp_vel(:, ip, 2) = this%prop_data%hamiltonian_elements(:, ip) * tmp_pos(:, ip, 1)
342 ! swap temporary variables
343 tmp_pos(:, ip, 1) = tmp_pos(:, ip, 2)
344 tmp_vel(:, ip, 1) = tmp_vel(:, ip, 2)
345 ! accumulate components of Taylor expansion
346 this%pos(:, ip) = this%pos(:, ip) + factor * tmp_pos(:, ip, 1)
347 this%vel(:, ip) = this%vel(:, ip) + factor * tmp_vel(:, ip, 1)
348 end do
349 end do
350 safe_deallocate_a(tmp_pos)
351 safe_deallocate_a(tmp_vel)
352 updated_quantities = ["position", "velocity"]
353
355 ! only correct for dt/2 if not converged yet
356 if (.not. this%is_tolerance_reached(prop%scf_tol)) then
357 this%pos(:, 1:this%np) = m_half*(this%pos(:, 1:this%np) + this%prop_data%save_pos(:, 1:this%np))
358 this%vel(:, 1:this%np) = m_half*(this%vel(:, 1:this%np) + this%prop_data%save_vel(:, 1:this%np))
359 updated_quantities = ["position", "velocity"]
360 end if
361
362 case default
363 done = .false.
364 end select
365
366 class default
367 done = .false.
368 end select
369
370 call profiling_out(trim(this%namespace%get())//":"//trim(operation%id))
373
374 ! ---------------------------------------------------------
375 logical function classical_particles_is_tolerance_reached(this, tol) result(converged)
376 class(classical_particles_t), intent(in) :: this
377 real(real64), intent(in) :: tol
378
379 integer :: ip
380 real(real64) :: change, max_change
381
383
384 ! Here we put the criterion that maximum acceleration change is below the tolerance
385 max_change = m_zero
386 do ip = 1, this%np
387 change = sum((this%prop_data%prev_tot_force(1:this%space%dim, ip) - this%tot_force(1:this%space%dim, ip))**2) / &
388 this%mass(ip)
389 if (change > max_change) then
390 max_change = change
391 end if
392 end do
393 converged = max_change < tol**2
394
395 write(message(1), '(a, e13.6, a, e13.6)') "Debug: -- Maximum change in acceleration ", &
396 sqrt(max_change), " and tolerance ", tol
397 call messages_info(1, namespace=this%namespace, debug_only=.true.)
398
401
402 ! ---------------------------------------------------------
403 subroutine classical_particles_update_quantity(this, label)
404 class(classical_particles_t), intent(inout) :: this
405 character(len=*), intent(in) :: label
406
408
409 select case (label)
410 case default
411 message(1) = "Incompatible quantity."
412 call messages_fatal(1, namespace=this%namespace)
413 end select
414
417
418 ! ---------------------------------------------------------
419 subroutine classical_particles_init_interaction_as_partner(partner, interaction)
420 class(classical_particles_t), intent(in) :: partner
421 class(interaction_surrogate_t), intent(inout) :: interaction
422
424
425 select type (interaction)
426 class default
427 message(1) = "Unsupported interaction."
428 call messages_fatal(1, namespace=partner%namespace)
429 end select
430
433
434 ! ---------------------------------------------------------
435 subroutine classical_particles_copy_quantities_to_interaction(partner, interaction)
436 class(classical_particles_t), intent(inout) :: partner
437 class(interaction_surrogate_t), intent(inout) :: interaction
438
440
441 select type (interaction)
442 class default
443 message(1) = "Unsupported interaction."
444 call messages_fatal(1, namespace=partner%namespace)
445 end select
446
449
450 ! ---------------------------------------------------------
452 class(classical_particles_t), intent(inout) :: this
453
455
456 ! Store previous force, as it is used as SCF criterium
457 select type (prop => this%algo)
458 class is (propagator_t)
459 if (prop%predictor_corrector) then
460 this%prop_data%prev_tot_force(1:this%space%dim, 1:this%np) = this%tot_force(1:this%space%dim, 1:this%np)
461 end if
462 end select
463
466
467 ! ---------------------------------------------------------
469 class(classical_particles_t), intent(inout) :: this
471 type(interaction_iterator_t) :: iter
472
474
475 ! Compute the total force acting on the classical particles
476 this%tot_force(1:this%space%dim, 1:this%np) = m_zero
477 call iter%start(this%interactions)
478 do while (iter%has_next())
479 select type (interaction => iter%get_next())
480 class is (force_interaction_t)
481 this%tot_force(1:this%space%dim, 1:this%np) = this%tot_force(1:this%space%dim, 1:this%np) + &
482 interaction%force(1:this%space%dim, 1:this%np)
483 end select
484 end do
485
488
489 ! ---------------------------------------------------------
490 subroutine classical_particles_output_start(this)
491 class(classical_particles_t), intent(inout) :: this
492
493 integer :: iteration
494 real(real64) :: dt
495
497
498 select type (algo => this%algo)
499 class is (propagator_t)
500 dt = algo%dt
501 end select
502
503 iteration = this%iteration%counter()
504 if (iteration > 0) then
505 ! add one if restarting because the output routine is only called at the end of the timestep
506 iteration = iteration + 1
507 end if
508 ! Create output handle
509 call io_mkdir('td.general', this%namespace)
510 if (this%grp%is_root()) then
511 call write_iter_init(this%output_handle(output_coordinates), iteration, dt, &
512 trim(io_workpath("td.general/coordinates", this%namespace)))
513 call write_iter_init(this%output_handle(output_energy), iteration, dt, &
514 trim(io_workpath("td.general/energy", this%namespace)))
515 end if
516
517 ! Output info for first iteration
518 if (iteration == 0) then
519 call this%output_write()
520 end if
521
524
525 ! ---------------------------------------------------------
526 subroutine classical_particles_output_finish(this)
527 class(classical_particles_t), intent(inout) :: this
528
531 if (this%grp%is_root()) then
532 call write_iter_end(this%output_handle(output_coordinates))
533 call write_iter_end(this%output_handle(output_energy))
534 end if
535
538
539 ! ---------------------------------------------------------
540 subroutine classical_particles_output_write(this)
541 class(classical_particles_t), intent(inout) :: this
542
543 integer :: idir, ii, ip
544 character(len=50) :: aux
545 real(real64) :: tmp(this%space%dim)
547 if (.not. this%grp%is_root()) return ! only first node outputs
548
550
551 if (this%iteration%counter() == 0) then
552 ! header
553 call write_iter_clear(this%output_handle(output_coordinates))
554 call write_iter_string(this%output_handle(output_coordinates), &
555 '############################################################################')
556 call write_iter_nl(this%output_handle(output_coordinates))
557 call write_iter_string(this%output_handle(output_coordinates),'# HEADER')
558 call write_iter_nl(this%output_handle(output_coordinates))
559
560 ! first line: column names
561 call write_iter_header_start(this%output_handle(output_coordinates))
562 do ip = 1, this%np
563 do idir = 1, this%space%dim
564 write(aux, '(a2,i3,a1,i3,a1)') 'x(', ip, ',', idir, ')'
565 call write_iter_header(this%output_handle(output_coordinates), aux)
566 end do
567 end do
568 do ip = 1, this%np
569 do idir = 1, this%space%dim
570 write(aux, '(a2,i3,a1,i3,a1)') 'v(', ip, ',', idir, ')'
571 call write_iter_header(this%output_handle(output_coordinates), aux)
572 end do
573 end do
574 do ip = 1, this%np
575 do idir = 1, this%space%dim
576 write(aux, '(a2,i3,a1,i3,a1)') 'f(', ip, ',', idir, ')'
577 call write_iter_header(this%output_handle(output_coordinates), aux)
578 end do
579 end do
580 call write_iter_nl(this%output_handle(output_coordinates))
581
582 ! second line: units
583 call write_iter_string(this%output_handle(output_coordinates), '#[Iter n.]')
584 call write_iter_header(this%output_handle(output_coordinates), '[' // trim(units_abbrev(units_out%time)) // ']')
585 do ip = 1, this%np
586 do idir = 1, this%space%dim
587 call write_iter_header(this%output_handle(output_coordinates), '[' // trim(units_abbrev(units_out%length)) // ']')
588 end do
589 end do
590 do ip = 1, this%np
591 do idir = 1, this%space%dim
592 call write_iter_header(this%output_handle(output_coordinates), '[' // trim(units_abbrev(units_out%velocity)) // ']')
593 end do
594 end do
595 do ip = 1, this%np
596 do idir = 1, this%space%dim
597 call write_iter_header(this%output_handle(output_coordinates), '[' // trim(units_abbrev(units_out%force)) // ']')
598 end do
599 end do
600 call write_iter_nl(this%output_handle(output_coordinates))
601
602 call write_iter_string(this%output_handle(output_coordinates), &
603 '############################################################################')
604 call write_iter_nl(this%output_handle(output_coordinates))
605 end if
606
607 call write_iter_start(this%output_handle(output_coordinates))
608
609 ! Position
610 do ip = 1, this%np
611 tmp(1:this%space%dim) = units_from_atomic(units_out%length, this%pos(1:this%space%dim, ip))
612 call write_iter_double(this%output_handle(output_coordinates), tmp, this%space%dim)
613 end do
614 ! Velocity
615 do ip = 1, this%np
616 tmp(1:this%space%dim) = units_from_atomic(units_out%velocity, this%vel(1:this%space%dim, ip))
617 call write_iter_double(this%output_handle(output_coordinates), tmp, this%space%dim)
618 end do
619 ! Force
620 do ip = 1, this%np
621 tmp(1:this%space%dim) = units_from_atomic(units_out%force, this%tot_force(1:this%space%dim, ip))
622 call write_iter_double(this%output_handle(output_coordinates), tmp, this%space%dim)
623 end do
624
625 call write_iter_nl(this%output_handle(output_coordinates))
626
627 ! Energies
628 if (this%iteration%counter() == 0) then
629 ! header
630 call write_iter_clear(this%output_handle(output_energy))
631 call write_iter_string(this%output_handle(output_energy), &
632 '############################################################################')
633 call write_iter_nl(this%output_handle(output_energy))
634 call write_iter_string(this%output_handle(output_energy),'# HEADER')
635 call write_iter_nl(this%output_handle(output_energy))
636
637 ! first line: column names
638 call write_iter_header_start(this%output_handle(output_energy))
639 call write_iter_header(this%output_handle(output_energy), 'Total')
640 call write_iter_header(this%output_handle(output_energy), 'Kinetic')
641 call write_iter_header(this%output_handle(output_energy), 'Potential')
642 call write_iter_header(this%output_handle(output_energy), 'Internal')
643 call write_iter_nl(this%output_handle(output_energy))
644
645 ! second line: units
646 call write_iter_string(this%output_handle(output_energy), '#[Iter n.]')
647 call write_iter_header(this%output_handle(output_energy), '[' // trim(units_abbrev(units_out%time)) // ']')
648 do ii = 1, 4
649 call write_iter_header(this%output_handle(output_energy), '[' // trim(units_abbrev(units_out%energy)) // ']')
650 end do
651 call write_iter_nl(this%output_handle(output_energy))
652
653 call write_iter_string(this%output_handle(output_energy), &
654 '############################################################################')
655 call write_iter_nl(this%output_handle(output_energy))
656 end if
657 call write_iter_start(this%output_handle(output_energy))
658
659 call write_iter_double(this%output_handle(output_energy), units_from_atomic(units_out%energy, this%total_energy), 1)
660 call write_iter_double(this%output_handle(output_energy), units_from_atomic(units_out%energy, this%kinetic_energy), 1)
661 call write_iter_double(this%output_handle(output_energy), units_from_atomic(units_out%energy, this%potential_energy), 1)
662 call write_iter_double(this%output_handle(output_energy), units_from_atomic(units_out%energy, this%internal_energy), 1)
663 call write_iter_nl(this%output_handle(output_energy))
664
667
668 ! ---------------------------------------------------------
670 class(classical_particles_t), intent(inout) :: this
671
672 type(restart_basic_t) :: restart
673 integer :: restart_file_unit
674 integer :: ierr
675
677
678 if (this%grp%is_root()) then
679 call write_iter_flush(this%output_handle(output_coordinates))
680 call write_iter_flush(this%output_handle(output_energy))
681 end if
682
683 call profiling_in(trim(this%namespace%get())//":"//"RESTART_WRITE")
684
685 call restart%init(this%namespace, restart_td, restart_type_dump, ierr)
686
687 restart_file_unit = restart%open("restart_classical_particles")
688 if (restart%do_i_write()) then
689
690 write(restart_file_unit, *) this%np
691 write(restart_file_unit, *) this%mass(:)
692 write(restart_file_unit, *) this%pos(:,:)
693 write(restart_file_unit, *) this%vel(:,:)
694 write(restart_file_unit, *) this%tot_force(:,:)
695
696 end if
697 call restart%close(restart_file_unit)
698
699
700 if (this%iteration%counter() > 0) then
701 ! only initialized after the first time step
702 call this%prop_data%restart_write(this%namespace, this%algo)
703 end if
704
705 message(1) = "Successfully wrote restart data for system "//trim(this%namespace%get())
706 call messages_info(1, namespace=this%namespace)
707
708 call restart%end()
709
710 call profiling_out(trim(this%namespace%get())//":"//"RESTART_WRITE")
713
714 ! ---------------------------------------------------------
715 logical function classical_particles_restart_read_data(this)
716 class(classical_particles_t), intent(inout) :: this
717
718 type(restart_basic_t) :: restart
719 integer :: restart_file_unit
720 integer :: ierr
721
723 call profiling_in(trim(this%namespace%get())//":"//"RESTART_READ")
724
725 call restart%init(this%namespace, restart_td, restart_type_load, ierr)
726
728 if (ierr==0) then
729 restart_file_unit = restart%open("restart_classical_particles")
730
731 if (restart_file_unit /= -1) then
732 read(restart_file_unit, *) this%np
733 read(restart_file_unit, *) this%mass(:)
734 read(restart_file_unit, *) this%pos(:,:)
735 read(restart_file_unit, *) this%vel(:,:)
736 read(restart_file_unit, *) this%tot_force(:,:)
737
738 call restart%close(restart_file_unit)
739 call this%prop_data%initialize(this%algo, this%space%dim, this%np)
740 if (this%iteration%counter() > 0) then
741 ! only initialized after the first time step
742 classical_particles_restart_read_data = this%prop_data%restart_read(this%namespace, this%algo)
743 else
745 end if
746 end if
747 end if
748 call restart%end()
749
751 message(1) = "Successfully read restart data for system "//trim(this%namespace%get())
752 ! namespace should be added here at some point
753 call messages_info(1)
754 else
755 message(1) = "Failed to read restart data for system "//trim(this%namespace%get())
756 ! namespace should be added here at some point
757 call messages_info(1)
758 end if
759
760 call profiling_out(trim(this%namespace%get())//":"//"RESTART_READ")
763
764 ! ---------------------------------------------------------
766 class(classical_particles_t), intent(inout) :: this
767
768 integer :: ip
769
771
772 this%kinetic_energy = m_zero
773 do ip = 1, this%np
774 this%kinetic_energy = this%kinetic_energy + m_half * this%mass(ip) * sum(this%vel(:, ip)**2)
775 end do
776
779
780 ! ---------------------------------------------------------
781 function classical_particles_center_of_mass(this, mask, pseudo) result(pos)
782 class(classical_particles_t), intent(in) :: this
783 logical, optional, intent(in) :: mask(:)
784 logical, optional, intent(in) :: pseudo
785 real(real64) :: pos(this%space%dim)
786
787 real(real64) :: mass, total_mass
788 integer :: ip
789
791
792 assert(.not. this%space%is_periodic())
793
794 pos = m_zero
795 total_mass = m_zero
796 mass = m_one
797 do ip = 1, this%np
798 if (present(mask)) then
799 if (.not. mask(ip)) cycle
800 end if
801 if (.not. optional_default(pseudo, .false.)) then
802 mass = this%mass(ip)
803 end if
804 pos = pos + mass*this%pos(:, ip)
805 total_mass = total_mass + mass
806 end do
807 pos = pos/total_mass
808
811
812 ! ---------------------------------------------------------
813 function classical_particles_center_of_mass_vel(this) result(vel)
814 class(classical_particles_t), intent(in) :: this
815 real(real64) :: vel(this%space%dim)
816
817 real(real64) :: mass, total_mass
818 integer :: ip
819
821
822 vel = m_zero
823 total_mass = m_zero
824 do ip = 1, this%np
825 mass = this%mass(ip)
826 total_mass = total_mass + mass
827 vel = vel + mass*this%vel(:, ip)
828 end do
829 vel = vel/total_mass
830
833
834 ! ---------------------------------------------------------
835 function classical_particles_center(this) result(pos)
836 class(classical_particles_t), intent(in) :: this
837 real(real64) :: pos(this%space%dim)
838
839 real(real64) :: xmin(this%space%dim), xmax(this%space%dim)
840 integer :: ip, idir
841
843
844 xmin = m_huge
845 xmax = -m_huge
846 do ip = 1, this%np
847 do idir = 1, this%space%dim
848 if (this%pos(idir, ip) > xmax(idir)) xmax(idir) = this%pos(idir, ip)
849 if (this%pos(idir, ip) < xmin(idir)) xmin(idir) = this%pos(idir, ip)
850 end do
851 end do
852
853 pos = (xmax + xmin)/m_two
854
856 end function classical_particles_center
857
858 ! ---------------------------------------------------------
859 subroutine classical_particles_axis_large(this, x, x2)
860 class(classical_particles_t), intent(in) :: this
861 real(real64), intent(out) :: x(this%space%dim)
862 real(real64), intent(out) :: x2(this%space%dim)
863
864 integer :: ip, jp
865 real(real64) :: rmax, r, r2
866
868
869 assert(.not. this%space%is_periodic())
870
871 ! first get the further apart atoms
872 rmax = -m_huge
873 do ip = 1, this%np
874 do jp = 1, this%np/2 + 1
875 r = norm2(this%pos(:, ip) - this%pos(:, jp))
876 if (r > rmax) then
877 rmax = r
878 x = this%pos(:, ip) - this%pos(:, jp)
879 end if
880 end do
881 end do
882 x = x /norm2(x)
883
884 ! now let us find out what is the second most important axis
885 rmax = -m_huge
886 do ip = 1, this%np
887 r2 = sum(x * this%pos(:, ip))
888 r = norm2(this%pos(:, ip) - r2*x)
889 if (r > rmax) then
890 rmax = r
891 x2 = this%pos(:, ip) - r2*x
892 end if
893 end do
894
896 end subroutine classical_particles_axis_large
897
898 ! ---------------------------------------------------------
901 subroutine classical_particles_axis_inertia(this, x, x2, pseudo)
902 class(classical_particles_t), intent(in) :: this
903 real(real64), intent(out) :: x(this%space%dim)
904 real(real64), intent(out) :: x2(this%space%dim)
905 logical, intent(in) :: pseudo
906
907 real(real64) :: mass, tinertia(this%space%dim, this%space%dim), eigenvalues(this%space%dim)
908 integer :: ii, jj, ip
909 type(unit_t) :: unit
910
912
913 assert(.not. this%space%is_periodic())
914
915 ! first calculate the inertia tensor
916 tinertia = m_zero
917 mass = m_one
918 do ip = 1, this%np
919 if (.not. pseudo) mass = this%mass(ip)
920 do ii = 1, this%space%dim
921 tinertia(ii, :) = tinertia(ii, :) - mass*this%pos(ii, ip)*this%pos(:, ip)
922 tinertia(ii, ii) = tinertia(ii, ii) + mass*sum(this%pos(:, ip)**2)
923 end do
924 end do
925
926 unit = units_out%length**2
927 ! note: we always use amu for atomic masses, so no unit conversion to/from atomic is needed.
928 if (pseudo) then
929 write(message(1),'(a)') 'Moment of pseudo-inertia tensor [' // trim(units_abbrev(unit)) // ']'
930 else
931 write(message(1),'(a)') 'Moment of inertia tensor [amu*' // trim(units_abbrev(unit)) // ']'
932 end if
933 call messages_info(1, namespace=this%namespace)
934 call output_tensor(tinertia, this%space%dim, unit, write_average = .true., namespace=this%namespace)
935
936 call lalg_eigensolve(this%space%dim, tinertia, eigenvalues)
937
938 write(message(1),'(a,6f25.6)') 'Eigenvalues: ', units_from_atomic(unit, eigenvalues)
939 call messages_info(1, namespace=this%namespace)
940
941 ! make a choice to fix the sign of the axis.
942 do ii = 1, 2
943 jj = maxloc(abs(tinertia(:,ii)), dim = 1)
944 if (tinertia(jj,ii) < m_zero) tinertia(:,ii) = -tinertia(:,ii)
945 end do
946 x = tinertia(:,1)
947 x2 = tinertia(:,2)
948
951
952 ! ---------------------------------------------------------
953 subroutine classical_particles_end(this)
954 class(classical_particles_t), intent(inout) :: this
955
957
958 safe_deallocate_a(this%mass)
959 safe_deallocate_a(this%pos)
960 safe_deallocate_a(this%vel)
961 safe_deallocate_a(this%tot_force)
962 safe_deallocate_a(this%fixed)
963 safe_deallocate_a(this%lj_epsilon)
964 safe_deallocate_a(this%lj_sigma)
965
966 call system_end(this)
967
969 end subroutine classical_particles_end
970
972
973!! Local Variables:
974!! mode: f90
975!! coding: utf-8
976!! End:
Writes to the corresponding file and adds one to the iteration. Must be called after write_iter_init(...
Definition: write_iter.F90:169
This module implements the basic elements defining algorithms.
Definition: algorithm.F90:143
subroutine classical_particles_output_finish(this)
subroutine classical_particles_axis_large(this, x, x2)
logical function classical_particles_is_tolerance_reached(this, tol)
subroutine classical_particles_output_write(this)
subroutine, public classical_particles_update_quantity(this, label)
logical function classical_particles_do_algorithmic_operation(this, operation, updated_quantities)
subroutine classical_particles_update_kinetic_energy(this)
subroutine classical_particles_axis_inertia(this, x, x2, pseudo)
This subroutine assumes that the origin of the coordinates is the center of mass of the system.
subroutine classical_particles_output_start(this)
logical function classical_particles_restart_read_data(this)
subroutine, public classical_particles_init_interaction_as_partner(partner, interaction)
subroutine, public classical_particles_end(this)
subroutine, public classical_particles_init(this, np)
The init routine is a module level procedure This has the advantage that different classes can have d...
subroutine, public classical_particles_copy(this, cp_in)
subroutine classical_particles_update_interactions_finish(this)
subroutine classical_particles_update_interactions_start(this)
subroutine, public classical_particles_copy_quantities_to_interaction(partner, interaction)
real(real64) function, dimension(this%space%dim) classical_particles_center_of_mass_vel(this)
subroutine, public classical_particles_init_interaction(this, interaction)
subroutine, public classical_particles_restart_write_data(this)
real(real64) function, dimension(this%space%dim) classical_particles_center_of_mass(this, mask, pseudo)
real(real64) function, dimension(this%space%dim) classical_particles_center(this)
real(real64), parameter, public m_two
Definition: global.F90:192
real(real64), parameter, public m_huge
Definition: global.F90:208
real(real64), parameter, public m_zero
Definition: global.F90:190
real(real64), parameter, public m_four
Definition: global.F90:194
real(real64), parameter, public m_epsilon
Definition: global.F90:206
real(real64), parameter, public m_half
Definition: global.F90:196
real(real64), parameter, public m_one
Definition: global.F90:191
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:317
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:360
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:416
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:600
This module implements the multisystem debug functionality.
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:625
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
character(len=algo_label_len), parameter, public beeman_start
character(len=algo_label_len), parameter, public beeman_finish
character(len=algo_label_len), parameter, public beeman_compute_acc
character(len=algo_label_len), parameter, public beeman_correct_pos
character(len=algo_label_len), parameter, public beeman_predict_pos
character(len=algo_label_len), parameter, public beeman_correct_vel
character(len=algo_label_len), parameter, public beeman_predict_vel
character(len=algo_label_len), parameter, public expmid_2step_predict_dt
character(len=algo_label_len), parameter, public expmid_2step_finish
character(len=algo_label_len), parameter, public expmid_2step_correct_dt_2
character(len=algo_label_len), parameter, public update_hamiltonian
character(len=algo_label_len), parameter, public expmid_2step_predict_dt_2
character(len=algo_label_len), parameter, public expmid_2step_start
This module implements the basic propagator framework.
Definition: propagator.F90:119
character(len=algo_label_len), parameter, public store_current_status
Definition: propagator.F90:171
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
integer, parameter, public restart_type_dump
Definition: restart.F90:183
integer, parameter, public restart_td
Definition: restart.F90:156
integer, parameter, public restart_type_load
Definition: restart.F90:183
This module implements the abstract system type.
Definition: system.F90:120
subroutine, public system_end(this)
Definition: system.F90:1149
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
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:120
subroutine, public output_tensor(tensor, ndim, unit, write_average, iunit, namespace)
Definition: utils.F90:245
Explicit interfaces to C functions, defined in write_iter_low.cc.
Definition: write_iter.F90:116
Descriptor of one algorithmic operation.
Definition: algorithm.F90:165
These class extend the list and list iterator to make an interaction list.
surrogate interaction class to avoid circular dependencies between modules.
Abstract class implementing propagators.
Definition: propagator.F90:140
Systems (system_t) can expose quantities that can be used to calculate interactions with other system...
Definition: quantity.F90:173
restart_basic_t stores the basic information about a restart object.
Definition: restart.F90:216
Abstract class for systems.
Definition: system.F90:174
int true(void)