41#include "../base/base_uses.f90"
48 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'wiener_process'
67 CHARACTER(LEN=40) :: name
68 INTEGER :: iparticle, iparticle_kind, &
69 iparticle_local, nparticle, &
70 nparticle_kind, nparticle_local
71 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :, :) :: seed
82 NULLIFY (work_section, force_env)
83 cpassert(
ASSOCIATED(md_env))
85 CALL get_md_env(md_env=md_env, force_env=force_env, para_env=para_env, &
91 CALL force_env_get(force_env, force_env_section=force_env_section, &
96 CALL cp_subsys_get(subsys=subsys, atomic_kinds=atomic_kinds, local_particles=local_particles, &
99 nparticle_kind = atomic_kinds%n_els
100 nparticle = particles%n_els
103 ALLOCATE (local_particles%local_particle_set(nparticle_kind))
105 DO iparticle_kind = 1, nparticle_kind
106 nparticle_local = local_particles%n_el(iparticle_kind)
107 ALLOCATE (local_particles%local_particle_set(iparticle_kind)%rng(nparticle_local))
108 DO iparticle_local = 1, nparticle_local
109 ALLOCATE (local_particles%local_particle_set(iparticle_kind)%rng(iparticle_local)%stream)
115 ALLOCATE (seed(3, 2, nparticle))
118 seed(:, :, 1) = subsys%seed(:, :)
120 DO iparticle = 2, nparticle
121 seed(:, :, iparticle) =
next_rng_seed(seed(:, :, iparticle - 1))
128 DO iparticle_kind = 1, nparticle_kind
129 nparticle_local = local_particles%n_el(iparticle_kind)
130 DO iparticle_local = 1, nparticle_local
131 iparticle = local_particles%list(iparticle_kind)%array(iparticle_local)
132 WRITE (unit=name, fmt=
"(A,I8)")
"Wiener process for particle", iparticle
134 local_particles%local_particle_set(iparticle_kind)%rng(iparticle_local)% &
136 extended_precision=.true., seed=seed(:, :, iparticle))
143 NULLIFY (work_section)
145 subsection_name=
"RNG_INIT")
146 CALL init_local_particle_set(distribution_1d=local_particles, &
147 nparticle_kind=nparticle_kind, &
148 work_section=work_section)
161 SUBROUTINE init_local_particle_set(distribution_1d, nparticle_kind, &
165 INTEGER,
INTENT(in) :: nparticle_kind
168 CHARACTER(LEN=rng_record_length) :: rng_record
169 INTEGER :: iparticle, iparticle_kind, &
170 iparticle_local, nparticle_local
175 cpassert(
ASSOCIATED(distribution_1d))
177 IF (
ASSOCIATED(work_section))
THEN
180 DO iparticle_kind = 1, nparticle_kind
181 nparticle_local = distribution_1d%n_el(iparticle_kind)
182 DO iparticle_local = 1, nparticle_local
183 iparticle = distribution_1d%list(iparticle_kind)%array(iparticle_local)
184 IF (iparticle == distribution_1d%list(iparticle_kind)%array(iparticle_local))
THEN
186 keyword_name=
"_DEFAULT_KEYWORD_", &
187 i_rep_val=iparticle, &
189 distribution_1d%local_particle_set(iparticle_kind)%rng(iparticle_local)% &
197 END SUBROUTINE init_local_particle_set
212 CHARACTER(LEN=40) :: name
214 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :, :) :: seed
215 REAL(kind=
dp),
DIMENSION(3, 2) :: initial_seed
217 IF (.NOT.
ASSOCIATED(meta_env))
RETURN
224 ALLOCATE (seed(3, 2, meta_env%n_colvar))
226 seed(:, :, 1) = initial_seed
227 DO i_c = 2, meta_env%n_colvar
235 DO i_c = 1, meta_env%n_colvar
236 WRITE (unit=name, fmt=
"(A,I8)")
"Wiener process for COLVAR", i_c
239 extended_precision=.true., seed=seed(:, :, i_c))
represent a simple array based list of the given type
types that represent a subsys, i.e. a part of the system
subroutine, public cp_subsys_get(subsys, ref_count, atomic_kinds, atomic_kind_set, particles, particle_set, local_particles, molecules, molecule_set, molecule_kinds, molecule_kind_set, local_molecules, para_env, colvar_p, shell_particles, core_particles, gci, multipoles, natom, nparticle, ncore, nshell, nkind, atprop, virial, results, cell)
returns information about various attributes of the given subsys
stores a lists of integer that are local to a processor. The idea is that these integers represent ob...
Interface for the force calculations.
recursive subroutine, public force_env_get(force_env, in_use, fist_env, qs_env, meta_env, fp_env, subsys, para_env, potential_energy, additional_potential, kinetic_energy, harmonic_shell, kinetic_shell, cell, sub_force_env, qmmm_env, qmmmx_env, eip_env, pwdft_env, globenv, input, force_env_section, method_name_id, root_section, mixed_env, nnp_env, embed_env, ipi_env)
returns various attributes about the force environment
Defines the basic variable types.
integer, parameter, public dp
pure logical function, public need_per_atom_wiener_process(md_env)
...
subroutine, public get_md_env(md_env, itimes, constant, used_time, cell, simpar, npt, force_env, para_env, reftraj, t, init, first_time, fe_env, thermostats, barostat, thermostat_coeff, thermostat_part, thermostat_shell, thermostat_baro, thermostat_fast, thermostat_slow, md_ener, averages, thermal_regions, ehrenfest_md)
get components of MD environment type
Interface to the message passing library MPI.
Parallel (pseudo)random number generator (RNG) for multiple streams and substreams of random numbers.
type(rng_stream_type) function, public rng_stream_type_from_record(rng_record)
Create a RNG stream from a record given as an internal file (string).
real(kind=dp) function, dimension(3, 2), public next_rng_seed(seed)
Get the seed for the next RNG stream w.r.t. a given seed.
integer, parameter, public rng_record_length
integer, parameter, public gaussian
represent a simple array based list of the given type
Type for storing MD parameters.
Utilities for string manipulations.
subroutine, public compress(string, full)
Eliminate multiple space characters in a string. If full is .TRUE., then all spaces are eliminated.
Handling of the Wiener process currently employed in turn of the Langevin dynamics.
subroutine, public create_wiener_process(md_env)
Create a Wiener process for Langevin dynamics and initialize an independent random number generator f...
subroutine, public create_wiener_process_cv(meta_env)
Create a Wiener process for Langevin dynamics used for metadynamics and initialize an independent ran...
represent a list of objects
represents a system: atoms, molecules, their pos,vel,...
structure to store local (to a processor) ordered lists of integers.
wrapper to abstract the force evaluation of the various methods
stores all the informations relevant to an mpi environment
represent a list of objects
Simulation parameter type for molecular dynamics.