Octopus
td_write_classical.F90
Go to the documentation of this file.
1! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22 use debug_oct_m
23 use global_oct_m
26 use io_oct_m
28 use mpi_oct_m
32 use space_oct_m
34 use unit_oct_m
37
38 implicit none
39
40 private
41
42 public :: &
45
46contains
47
48 ! ---------------------------------------------------------
49 subroutine td_write_coordinates(out_coords, natoms, space, pos, vel, tot_forces, iter)
50 type(c_ptr), intent(inout) :: out_coords
51 integer, intent(in) :: natoms
52 class(space_t), intent(in) :: space
53 real(real64), intent(in) :: pos(:,:)
54 real(real64), intent(in) :: vel(:,:)
55 real(real64), intent(in) :: tot_forces(:,:)
56 integer, intent(in) :: iter
57
58 integer :: iatom, idir
59 character(len=50) :: aux
60 real(real64) :: tmp(space%dim)
61
62 if (.not. mpi_grp_is_root(mpi_world)) return ! only first node outputs
63
64 push_sub(td_write_coordinates)
65
66 if (iter == 0) then
67 call td_write_print_header_init(out_coords)
68
69 ! first line: column names
70 call write_iter_header_start(out_coords)
71
72 do iatom = 1, natoms
73 do idir = 1, space%dim
74 write(aux, '(a2,i0,a1,i1,a1)') 'x(', iatom, ',', idir, ')'
75 call write_iter_header(out_coords, aux)
76 end do
77 end do
78 do iatom = 1, natoms
79 do idir = 1, space%dim
80 write(aux, '(a2,i0,a1,i1,a1)') 'v(', iatom, ',', idir,')'
81 call write_iter_header(out_coords, aux)
82 end do
83 end do
84 do iatom = 1, natoms
85 do idir = 1, space%dim
86 write(aux, '(a2,i0,a1,i1,a1)') 'f(', iatom, ',', idir,')'
87 call write_iter_header(out_coords, aux)
88 end do
89 end do
90 call write_iter_nl(out_coords)
91
92 ! second line: units
93 call write_iter_string(out_coords, '#[Iter n.]')
94 call write_iter_header(out_coords, '[' // trim(units_abbrev(units_out%time)) // ']')
95 call write_iter_string(out_coords, &
96 'Positions in ' // trim(units_abbrev(units_out%length)) // &
97 ', Velocities in '// trim(units_abbrev(units_out%velocity)) // &
98 ', Forces in ' // trim(units_abbrev(units_out%force)))
99 call write_iter_nl(out_coords)
100
101 call td_write_print_header_end(out_coords)
102 end if
103
104 call write_iter_start(out_coords)
105
106 do iatom = 1, natoms
107 tmp(1:space%dim) = units_from_atomic(units_out%length, pos(:, iatom))
108 call write_iter_double(out_coords, tmp, space%dim)
109 end do
110 do iatom = 1, natoms
111 tmp(1:space%dim) = units_from_atomic(units_out%velocity, vel(:, iatom))
112 call write_iter_double(out_coords, tmp, space%dim)
113 end do
114 do iatom = 1, natoms
115 tmp(1:space%dim) = units_from_atomic(units_out%force, tot_forces(:, iatom))
116 call write_iter_double(out_coords, tmp, space%dim)
117 end do
118 call write_iter_nl(out_coords)
119
120 pop_sub(td_write_coordinates)
121 end subroutine td_write_coordinates
122
123 ! ---------------------------------------------------------
124 subroutine td_write_sep_coordinates(out_coords, natoms, space, pos, vel, tot_forces, iter, which)
125 type(c_ptr), intent(inout) :: out_coords
126 integer, intent(in) :: natoms
127 class(space_t), intent(in) :: space
128 real(real64), intent(in) :: pos(:,:)
129 real(real64), intent(in) :: vel(:,:)
130 real(real64), intent(in) :: tot_forces(:,:)
131 integer, intent(in) :: iter
132 integer, intent(in) :: which !1=xyz, 2=velocity, 3=force
133
134 integer, parameter :: COORDINATES=1
135 integer, parameter :: VELOCITIES=2
136 integer, parameter :: FORCES=3
137 integer :: iatom, idir
138 character(len=50) :: aux
139 real(real64) :: tmp(space%dim)
140
141 if (.not. mpi_grp_is_root(mpi_world)) return ! only first node outputs
144
145 if (iter == 0) then
146 call td_write_print_header_init(out_coords)
147
148 ! first line: column names
149 call write_iter_header_start(out_coords)
150
151 do iatom = 1, natoms
152 do idir = 1, space%dim
153 select case (which)
154 case (coordinates)
155 write(aux, '(a2,i0,a1,i1,a1)') 'x(', iatom, ',', idir, ')'
156 case (velocities)
157 write(aux, '(a2,i0,a1,i1,a1)') 'v(', iatom, ',', idir,')'
158 case (forces)
159 write(aux, '(a2,i0,a1,i1,a1)') 'f(', iatom, ',', idir,')'
160 case default
161 assert(.false.)
162 end select
163 call write_iter_header(out_coords, aux)
164 end do
165 end do
166 call write_iter_nl(out_coords)
167
168 ! second line: units
169 call write_iter_string(out_coords, '#[Iter n.]')
170 call write_iter_header(out_coords, '[' // trim(units_abbrev(units_out%time)) // ']')
171 select case (which)
172 case (coordinates)
173 call write_iter_string(out_coords, &
174 'Positions in ' // trim(units_abbrev(units_out%length)))
175 case (velocities)
176 call write_iter_string(out_coords, &
177 'Velocities in ' // trim(units_abbrev(units_out%velocity)))
178 case (forces)
179 call write_iter_string(out_coords, &
180 'Forces in ' // trim(units_abbrev(units_out%force)))
181 end select
182 call write_iter_nl(out_coords)
183
184 call td_write_print_header_end(out_coords)
185 end if
186
187 call write_iter_start(out_coords)
188
189 select case (which)
190 case (coordinates)
191 do iatom = 1, natoms
192 tmp(1:space%dim) = units_from_atomic(units_out%length, pos(:, iatom))
193 call write_iter_double(out_coords, tmp, space%dim)
194 end do
195 case (velocities)
196 do iatom = 1, natoms
197 tmp(1:space%dim) = units_from_atomic(units_out%velocity, vel(:, iatom))
198 call write_iter_double(out_coords, tmp, space%dim)
199 end do
200 case (forces)
201 do iatom = 1, natoms
202 tmp(1:space%dim) = units_from_atomic(units_out%force, tot_forces(:, iatom))
203 call write_iter_double(out_coords, tmp, space%dim)
204 end do
205 end select
206
207 call write_iter_nl(out_coords)
208
210 end subroutine td_write_sep_coordinates
211
213
214!! Local Variables:
215!! mode: f90
216!! coding: utf-8
217!! End:
218
Writes to the corresponding file and adds one to the iteration. Must be called after write_iter_init(...
Definition: write_iter.F90:167
Definition: io.F90:114
This module contains some common usage patterns of MPI routines.
Definition: mpi_lib.F90:115
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
this module contains the low-level part of the output system
Definition: output_low.F90:115
subroutine, public td_write_coordinates(out_coords, natoms, space, pos, vel, tot_forces, iter)
subroutine, public td_write_sep_coordinates(out_coords, natoms, space, pos, vel, tot_forces, iter, which)
subroutine, public td_write_print_header_init(out)
subroutine, public td_write_print_header_end(out)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
Explicit interfaces to C functions, defined in write_iter_low.cc.
Definition: write_iter.F90:114