Octopus
classical_particle.F90
Go to the documentation of this file.
1!! Copyright (C) 2019 N. Tancogne-Dejean
2!! Copyright (C) 2020 M. Oliveira
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
25 use debug_oct_m
27 use global_oct_m
32 use io_oct_m
35 use mpi_oct_m
37 use parser_oct_m
42 use space_oct_m
43 use system_oct_m
44 use unit_oct_m
46
47 implicit none
48
49 private
50 public :: &
53
54 integer, parameter :: &
55 OUTPUT_COORDINATES = 1, &
57
61 contains
62 procedure :: init_interaction => classical_particle_init_interaction
63 procedure :: initialize => classical_particle_initialize
64 procedure :: update_quantity => classical_particle_update_quantity
65 procedure :: init_interaction_as_partner => classical_particle_init_interaction_as_partner
66 procedure :: copy_quantities_to_interaction => classical_particle_copy_quantities_to_interaction
69
70 interface classical_particle_t
71 procedure classical_particle_constructor
72 end interface classical_particle_t
73
74contains
75
76 ! ---------------------------------------------------------
81 function classical_particle_constructor(namespace, grp) result(sys)
82 class(classical_particle_t), pointer :: sys
83 type(namespace_t), intent(in) :: namespace
84 type(mpi_grp_t), intent(in) :: grp
85
87
88 allocate(sys)
89
90 call classical_particle_init(sys, namespace, grp)
91
94
95 ! ---------------------------------------------------------
100 ! ---------------------------------------------------------
101 subroutine classical_particle_init(this, namespace, grp)
102 class(classical_particle_t), intent(inout) :: this
103 type(namespace_t), intent(in) :: namespace
104 type(mpi_grp_t), intent(in) :: grp
105
107
108 this%namespace = namespace
109 this%space = space_t(namespace)
110 call this%init_parallelization(grp)
111
112 if (this%space%is_periodic()) then
113 call messages_not_implemented('Classical particle for periodic systems')
114 end if
115
116 call messages_print_with_emphasis(msg="Classical Particle", namespace=namespace)
118 call classical_particles_init(this, 1)
119
120 !%Variable ParticleMass
121 !%Type float
122 !%Section ClassicalParticles
123 !%Description
124 !% Mass of classical particle in Kg.
125 !%End
126 call parse_variable(namespace, 'ParticleMass', m_one, this%mass(1))
127 call messages_print_var_value('ParticleMass', this%mass(1), namespace=namespace)
128
129 this%supported_interactions = [this%supported_interactions, gravity, lennard_jones]
130 this%supported_interactions_as_partner = [this%supported_interactions_as_partner, gravity, lennard_jones]
131
132 !%Variable LennardJonesEpsilon
133 !%Type float
134 !%Section ClassicalParticles
135 !%Description
136 !% Epsilon parameter (dispersion energy) of Lennard-Jones interaction for this species.
137 !% In case two particles have a different epsilon, the combination rule will be computed
138 !% <math>\epsilon_{12} = \sqrt{\epsilon_1 + \epsilon_2}</math>.
139 !%End
140 call parse_variable(namespace, 'LennardJonesEpsilon', m_one, this%lj_epsilon(1))
141 call messages_print_var_value('LennardJonesEpsilon', this%lj_epsilon(1), namespace=namespace)
142
143 !%Variable LennardJonesSigma
144 !%Type float
145 !%Section ClassicalParticles
146 !%Description
147 !% Sigma parameter (particle size) of Lennard-Jones interaction for this species.
148 !% In case two particles have a different sigma, the combination rule will be computed
149 !% <math>\sigma_{12} = (\sigma_1 + \sigma_2) / 2</math>.
150 !%End
151 call parse_variable(namespace, 'LennardJonesSigma', m_one, this%lj_sigma(1))
152 call messages_print_var_value('LennardJonesSigma', this%lj_sigma(1), namespace=namespace)
153
154 call messages_print_with_emphasis(namespace=namespace)
159 ! ---------------------------------------------------------
160 subroutine classical_particle_init_interaction(this, interaction)
161 class(classical_particle_t), target, intent(inout) :: this
162 class(interaction_t), intent(inout) :: interaction
163
165
166 select type (interaction)
167 type is (gravity_t)
168 call interaction%init(this%space%dim, 1, this%mass, this%pos)
169 type is (lennard_jones_t)
170 if (.not. (parse_is_defined(this%namespace, 'LennardJonesSigma') .and. &
171 parse_is_defined(this%namespace, 'LennardJonesEpsilon') )) then
172 write(message(1),'(a,es9.2)') 'Using default value for Lennard-Jones parameter.'
173 call messages_warning(1, namespace=this%namespace)
174 end if
175
176 call interaction%init(this%space%dim, 1, this%pos, this%lj_epsilon(1), this%lj_sigma(1))
177 class default
178 call classical_particles_init_interaction(this, interaction)
179 end select
180
183
184 ! ---------------------------------------------------------
185 subroutine classical_particle_initialize(this)
186 class(classical_particle_t), intent(inout) :: this
187
188 integer :: n_rows, idir
189 type(block_t) :: blk
190
191 real(real64) :: random_pos(1:this%space%dim), width_pos(1:this%space%dim)
192 real(real64) :: random_vel(1:this%space%dim), width_vel(1:this%space%dim)
193 integer(int64) :: seed
194 logical :: in_ensemble
195
197
198 in_ensemble = parse_is_defined(this%namespace, "NumberOfReplicas")
199
200 !%Variable ParticleInitialPosition
201 !%Type block
202 !%Section ClassicalParticles
203 !%Description
204 !% Initial position of classical particle, in Km.
205 !%End
206 this%pos = m_zero
207 if (parse_block(this%namespace, 'ParticleInitialPosition', blk) == 0) then
208 n_rows = parse_block_n(blk)
209 if (n_rows > 1) call messages_input_error(this%namespace, 'ParticleInitialPosition')
210
211 do idir = 1, this%space%dim
212 call parse_block_float(blk, 0, idir - 1, this%pos(idir, 1))
213 end do
214 call parse_block_end(blk)
215 end if
216 call messages_print_var_value('ParticleInitialPosition', this%pos(1:this%space%dim, 1), namespace=this%namespace)
217
218 !%Variable ParticleInitialVelocity
219 !%Type block
220 !%Section ClassicalParticles
221 !%Description
222 !% Initial velocity of classical particle in Km/s.
223 !%End
224 this%vel = m_zero
225 if (parse_block(this%namespace, 'ParticleInitialVelocity', blk) == 0) then
226 n_rows = parse_block_n(blk)
227 if (n_rows > 1) call messages_input_error(this%namespace, 'ParticleInitialVelocity')
228 do idir = 1, this%space%dim
229 call parse_block_float(blk, 0, idir - 1, this%vel(idir, 1))
230 end do
231 call parse_block_end(blk)
232 end if
233 call messages_print_var_value('ParticleInitialVelocity', this%vel(1:this%space%dim, 1), namespace=this%namespace)
234
235 if (in_ensemble) then
236
237 seed = this%namespace%get_hash32()
238 width_pos = m_zero
239 width_vel = m_zero
240
241 !%Variable ParticlePositionDistributionWidth
242 !%Type block
243 !%Section ClassicalParticles
244 !%Description
245 !% Width of the position random displacements of classical particle in Km.
246 !%End
247 if (parse_block(this%namespace, 'ParticlePositionDistributionWidth', blk) == 0) then
248 n_rows = parse_block_n(blk)
249 if (n_rows > 1) call messages_input_error(this%namespace, 'ParticlePositionDistributionWidth')
250 do idir = 1, this%space%dim
251 call parse_block_float(blk, 0, idir - 1, width_pos(idir))
252 end do
253 call parse_block_end(blk)
254 end if
255 call messages_print_var_value('ParticlePositionDistributionWidth', width_pos(1:this%space%dim), namespace=this%namespace)
256
257 !%Variable ParticleVelocityDistributionWidth
258 !%Type block
259 !%Section ClassicalParticles
260 !%Description
261 !% Width of the velocity random distribution of classical particle in Km/s.
262 !%End
263 if (parse_block(this%namespace, 'ParticleVelocityDistributionWidth', blk) == 0) then
264 n_rows = parse_block_n(blk)
265 if (n_rows > 1) call messages_input_error(this%namespace, 'ParticleVelocityDistributionWidth')
266 do idir = 1, this%space%dim
267 call parse_block_float(blk, 0, idir - 1, width_vel(idir))
268 end do
269 call parse_block_end(blk)
270 end if
271 call messages_print_var_value('ParticleVelocityDistributionWidth', width_vel(1:this%space%dim), namespace=this%namespace)
272
273 call shiftseed(seed, 1000_int64)
274 call quickrnd(seed, this%space%dim, random_pos)
275 call quickrnd(seed, this%space%dim, random_vel)
276
277
278 this%pos(:,1) = this%pos(:,1) + (random_pos - 0.5_real64) * width_pos
279 this%vel(:,1) = this%vel(:,1) + (random_vel - 0.5_real64) * width_vel
281 end if
282
284 end subroutine classical_particle_initialize
285
286 ! ---------------------------------------------------------
287 subroutine classical_particle_update_quantity(this, label)
288 class(classical_particle_t), intent(inout) :: this
289 character(len=*), intent(in) :: label
290
292
293 select case (label)
294 case default
295 ! Other quantities should be handled by the parent class
297 end select
298
301
302 ! ---------------------------------------------------------
303 subroutine classical_particle_init_interaction_as_partner(partner, interaction)
304 class(classical_particle_t), intent(in) :: partner
305 class(interaction_surrogate_t), intent(inout) :: interaction
306
308
309 select type (interaction)
310 type is (gravity_t)
311 interaction%partner_np = 1
312 safe_allocate(interaction%partner_mass(1))
313 safe_allocate(interaction%partner_pos(1:partner%space%dim, 1))
314 interaction%partner_pos = m_zero
315
316 type is (lennard_jones_t)
317 interaction%partner_np = 1
318 safe_allocate(interaction%partner_pos(1:partner%space%dim, 1))
319 ! in case the LennardJones epsilon and sigma of system and partner are different,
320 ! we compute the combination rules with arithmetic and geometric means:
321 ! (they give back the original parameters if they happen to be equal):
322 ! sigma_12 = (sigma_1 + sigma_2)/2, epsilon_12 = sqrt(epsilon_1 * epsilon_2)
323 interaction%lj_sigma = m_half * (partner%lj_sigma(1) + interaction%lj_sigma)
324 interaction%lj_epsilon = sqrt(partner%lj_epsilon(1) * interaction%lj_epsilon)
325
326 class default
327 ! Other interactions should be handled by the parent class
328 call classical_particles_init_interaction_as_partner(partner, interaction)
329 end select
330
333
334 ! ---------------------------------------------------------
335 subroutine classical_particle_copy_quantities_to_interaction(partner, interaction)
336 class(classical_particle_t), intent(inout) :: partner
337 class(interaction_surrogate_t), intent(inout) :: interaction
338
340
341 select type (interaction)
342 type is (gravity_t)
343 interaction%partner_mass(1) = partner%mass(1)
344 interaction%partner_pos(:,1) = partner%pos(:,1)
345
346 type is (lennard_jones_t)
347 interaction%partner_pos(:,1) = partner%pos(:,1)
348
349 class default
350 message(1) = "Unsupported interaction."
351 call messages_fatal(1, namespace=partner%namespace)
352 end select
353
356
357 ! ---------------------------------------------------------
358 subroutine classical_particle_finalize(this)
359 type(classical_particle_t), intent(inout) :: this
360
362
363 call classical_particles_end(this)
364
366 end subroutine classical_particle_finalize
367
369
370!! Local Variables:
371!! mode: f90
372!! coding: utf-8
373!! End:
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:182
This module implements the basic elements defining algorithms.
Definition: algorithm.F90:143
subroutine classical_particle_init_interaction_as_partner(partner, interaction)
subroutine, public classical_particle_init(this, namespace, grp)
The init routine is a module level procedure This has the advantage that different classes can have d...
class(classical_particle_t) function, pointer classical_particle_constructor(namespace, grp)
The factory routine (or constructor) allocates a pointer of the corresponding type and then calls the...
subroutine classical_particle_finalize(this)
integer, parameter output_energy
subroutine classical_particle_initialize(this)
subroutine classical_particle_update_quantity(this, label)
subroutine classical_particle_init_interaction(this, interaction)
subroutine classical_particle_copy_quantities_to_interaction(partner, interaction)
subroutine, public classical_particles_update_quantity(this, label)
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_init_interaction(this, interaction)
This module provides a random number generator for a normalized gaussian distribution.
real(real64), parameter, public m_zero
Definition: global.F90:200
real(real64), parameter, public m_half
Definition: global.F90:206
real(real64), parameter, public m_one
Definition: global.F90:201
integer, parameter, public lennard_jones
integer, parameter, public gravity
This module defines the abstract interaction_t class, and some auxiliary classes for interactions.
Definition: io.F90:116
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:898
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1068
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
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:463
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:623
This module implements the basic propagator framework.
Definition: propagator.F90:119
This module defines the quantity_t class and the IDs for quantities, which can be exposed by a system...
Definition: quantity.F90:140
subroutine, public shiftseed(iseed, n)
Definition: quickrnd.F90:260
This module implements the abstract system type.
Definition: system.F90:120
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
This module defines the unit system, used for input and output.
class for a neutral classical particle
Gravity interaction between two systems of particles. This should be used for testing purposes only....
Definition: gravity.F90:136
abstract interaction class
surrogate interaction class to avoid circular dependencies between modules.
Lennard-Jones interaction between two systems of particles.