Octopus
centergeom.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
21program centergeom
23 use global_oct_m
24 use io_oct_m
25 use ions_oct_m
27 use mpi_oct_m
29 use parser_oct_m
31 use unit_oct_m
33
34 implicit none
35
36 integer :: ierr
37
39
40 call getopt_init(ierr)
41 if (ierr == 0) call getopt_center_geom()
42 call getopt_end()
43
44 call parser_init()
45
46 call messages_init()
47
48 call io_init()
51
52 call center_geo()
53
55 call io_end()
56 call messages_end()
57
58 call parser_end()
59
60 call global_end()
61
62contains
63
64 subroutine center_geo()
65
66 integer, parameter :: &
67 NONE = 0, &
68 inertia = 1, &
69 pseudo = 2, &
70 large = 3
71
72 type(ions_t), pointer :: ions
73 real(real64), allocatable :: center(:), x1(:), x2(:), to(:)
74 integer :: axis_type, idir, default
75 type(block_t) :: blk
76
78
79 if (ions%space%is_periodic()) then
80 message(1) = "oct-center-geom utility does not work for periodic systems."
81 call messages_fatal(1, namespace=ions%namespace)
82 end if
83
84 safe_allocate(center(1:ions%space%dim))
85 safe_allocate(x1(1:ions%space%dim))
86 safe_allocate(x2(1:ions%space%dim))
87 safe_allocate(to(1:ions%space%dim))
88
89 if (parse_is_defined(ions%namespace, 'BoxCenter')) then
90 call messages_not_implemented("The centergeom utitlity assumes the center of the box "// &
91 "at the origin and is not compatible with a shifted BoxCenter")
92 end if
93
94 ! is there something to do
95 if (ions%natoms > 1) then
96
97 !%Variable MainAxis
98 !%Type block
99 !%Section Utilities::oct-center-geom
100 !%Description
101 !% A vector of reals defining the axis to which the molecule
102 !% should be aligned. If not present, the default value will
103 !% be the x-axis. For example in 3D:
104 !% <tt>
105 !% <br>%MainAxis
106 !% <br> 1 | 0 | 0
107 !% <br>%</tt>
108 !%End
109 if (parse_block(ions%namespace, 'MainAxis', blk) == 0) then
110 do idir = 1, ions%space%dim
111 call parse_block_float(blk, 0, idir - 1, to(idir))
112 end do
113 call parse_block_end(blk)
114 else
115 to(:) = m_zero
116 to(1) = m_one
117 end if
118 if (norm2(to) > m_epsilon) to = to / norm2(to)
119
120 write(message(1),'(a,6f15.6)') 'Using main axis ', to
121 call messages_info(1)
122
123 !%Variable AxisType
124 !%Type integer
125 !%Default inertia
126 !%Section Utilities::oct-center-geom
127 !%Description
128 !% After the structure is centered, it is also aligned to a set of orthogonal axes.
129 !% This variable decides which set of axes to use. Only implemented for 3D, in which case
130 !% the default is <tt>inertia</tt>; otherwise <tt>none</tt> is default and the only legal value.
131 !%Option none 0
132 !% Do not rotate. Will still give output regarding center of mass and moment of inertia.
133 !%Option inertia 1
134 !% The axis of inertia.
135 !%Option pseudo_inertia 2
136 !% Pseudo-axis of inertia, calculated considering all species to have equal mass.
137 !%Option large_axis 3
138 !% The larger axis of the molecule.
139 !%End
140 if (ions%space%dim == 3) then
141 default = inertia
142 else
143 default = none
144 end if
145 call parse_variable(ions%namespace, 'AxisType', default, axis_type)
146 call messages_print_var_option("AxisType", axis_type, namespace=ions%namespace)
147
148 if (ions%space%dim /= 3 .and. axis_type /= none) then
149 call messages_not_implemented("alignment to axes (AxisType /= none) in other than 3 dimensions", namespace=ions%namespace)
150 end if
151
152 select case (axis_type)
153 case (none, inertia, pseudo)
154 center = ions%center_of_mass(pseudo = (axis_type == pseudo))
155
156 write(message(1),'(3a,99f15.6)') 'Center of mass [', trim(units_abbrev(units_out%length)), '] = ', &
157 units_from_atomic(units_out%length, center)
158 call messages_info(1)
160 call ions%translate(center)
161 call ions%axis_inertia(x1, x2, pseudo = (axis_type == pseudo))
162 case (large)
163 center = ions%center()
164
165 write(message(1),'(3a,99f15.6)') 'Center [', trim(units_abbrev(units_out%length)), '] = ', &
166 units_from_atomic(units_out%length, center)
167 call messages_info(1)
168
169 call ions%translate(center)
170 call ions%axis_large(x1, x2)
171 case default
172 write(message(1), '(a,i2,a)') 'AxisType = ', axis_type, ' not known by Octopus.'
173 call messages_fatal(1, namespace=ions%namespace)
174 end select
175
176 write(message(1),'(a,99f15.6)') 'Found primary axis ', x1
177 write(message(2),'(a,99f15.6)') 'Found secondary axis ', x2
178 call messages_info(2)
179
180 if (axis_type /= none) then
181 call ions%rotate(x1, x2, to)
182 end if
183
184 end if
185
186 ! recenter
187 center = ions%center()
188 call ions%translate(center)
189
190 ! write adjusted geometry
191 call ions%write_xyz('./adjusted')
192
193 safe_deallocate_a(center)
194 safe_deallocate_a(x1)
195 safe_deallocate_a(x2)
196 safe_deallocate_a(to)
197 safe_deallocate_p(ions)
198
199 end subroutine center_geo
200
201end program centergeom
202
203!! Local Variables:
204!! mode: f90
205!! coding: utf-8
206!! End:
subroutine center_geo()
Definition: centergeom.F90:160
program centergeom
Definition: centergeom.F90:116
subroutine, public getopt_init(ierr)
Initializes the getopt machinery. Must be called before attempting to parse the options....
subroutine, public getopt_end
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
Definition: global.F90:494
real(real64), parameter, public m_zero
Definition: global.F90:200
type(mpi_comm), parameter, public serial_dummy_comm
Alias MPI_COMM_UNDEFINED for the specific use case of initialising Octopus utilities with no MPI supp...
Definition: global.F90:294
subroutine, public init_octopus_globals(comm)
Initialise Octopus-specific global constants and files. This routine performs no initialisation calls...
Definition: global.F90:432
real(real64), parameter, public m_epsilon
Definition: global.F90:216
real(real64), parameter, public m_one
Definition: global.F90:201
Definition: io.F90:116
subroutine, public io_init(defaults)
If the argument defaults is present and set to true, then the routine will not try to read anything f...
Definition: io.F90:165
subroutine, public io_end()
Definition: io.F90:271
subroutine, public messages_end()
Definition: messages.F90:273
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1068
subroutine, public messages_init(output_dir)
Definition: messages.F90:220
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_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
type(namespace_t), public global_namespace
Definition: namespace.F90:135
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:463
subroutine, public parser_init()
Initialise the Octopus parser.
Definition: parser.F90:410
subroutine, public parser_end()
End the Octopus parser.
Definition: parser.F90:442
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:623
subroutine, public profiling_end(namespace)
Definition: profiling.F90:415
subroutine, public profiling_init(namespace)
Create profiling subdirectory.
Definition: profiling.F90:257
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
subroutine, public unit_system_init(namespace)