62#include "./base/base_uses.f90"
68 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'input_restart_force_eval'
85 write_binary_restart_file, respa)
89 LOGICAL,
INTENT(IN) :: write_binary_restart_file
90 LOGICAL,
OPTIONAL :: respa
92 INTEGER :: i, i_subsys, iforce_eval, myid, &
94 INTEGER,
DIMENSION(:),
POINTER :: i_force_eval
95 LOGICAL :: is_present, multiple_subsys, &
97 REAL(kind=
dp),
DIMENSION(:),
POINTER :: work
102 force_env_sections, qmmm_section, &
103 rng_section, subsys_section, &
107 NULLIFY (rng_section, subsys_section, cell_section, virial, subsys, cell, dft_control, work, tmp_section)
117 CALL multiple_fe_list(force_env_sections, root_section, i_force_eval, nforce_eval)
122 i_rep_section=i_force_eval(iforce_eval))
123 CALL update_subsys(subsys_section, force_env, skip_vel_section, &
124 write_binary_restart_file)
127 IF (
PRESENT(respa))
THEN
130 iforce_eval = i_subsys
132 i_rep_section=i_force_eval(iforce_eval))
133 CALL update_subsys(subsys_section, force_env, skip_vel_section, &
134 write_binary_restart_file)
140 CALL update_rng_particle(rng_section, force_env)
143 i_rep_section=i_force_eval(iforce_eval))
144 CALL update_qmmm(qmmm_section, force_env)
147 IF (nforce_eval == 1)
THEN
148 IF (
ASSOCIATED(force_env%qs_env))
THEN
149 CALL get_qs_env(force_env%qs_env, dft_control=dft_control)
150 IF (dft_control%qs_control%cdft)
THEN
152 i_rep_section=i_force_eval(iforce_eval))
154 i_val=dft_control%qs_control%cdft_control%ienergy)
155 ALLOCATE (work(
SIZE(dft_control%qs_control%cdft_control%strength)))
156 work = dft_control%qs_control%cdft_control%strength
164 IF (nforce_eval > 1)
THEN
167 CALL section_vals_val_get(root_section,
"MULTIPLE_FORCE_EVALS%MULTIPLE_SUBSYS", l_val=multiple_subsys)
168 IF (virial%pv_availability .AND. multiple_subsys)
THEN
169 DO iforce_eval = 2, nforce_eval
171 i_rep_section=i_force_eval(iforce_eval))
173 CALL update_cell_section(cell, cell_section)
177 IF (
ASSOCIATED(force_env%mixed_env))
THEN
178 IF (force_env%mixed_env%do_mixed_cdft)
THEN
179 DO iforce_eval = 2, nforce_eval
181 i_rep_section=i_force_eval(iforce_eval))
182 ALLOCATE (work(
SIZE(force_env%mixed_env%strength, 2)))
183 work = force_env%mixed_env%strength(iforce_eval - 1, :)
187 i_val=force_env%mixed_env%cdft_control%sim_step)
193 IF (
ASSOCIATED(force_env%qs_env))
THEN
194 CALL get_qs_env(force_env%qs_env, dft_control=dft_control)
196 i_rep_section=i_force_eval(iforce_eval))
202 IF (
ASSOCIATED(dft_control%efield_fields))
THEN
205 DO i = 1,
SIZE(dft_control%efield_fields)
206 work = dft_control%efield_fields(1)%efield%vec_pot_initial
213 DEALLOCATE (i_force_eval)
225 SUBROUTINE update_qmmm(qmmm_section, force_env)
231 REAL(kind=
dp),
DIMENSION(:),
POINTER :: work
233 SELECT CASE (force_env%in_use)
239 work = force_env%qmmm_env%qm%transl_v
243 END SUBROUTINE update_qmmm
254 SUBROUTINE update_rng_particle(rng_section, force_env)
259 CHARACTER(LEN=rng_record_length) :: rng_record
260 INTEGER :: iparticle, iparticle_kind, &
261 iparticle_local, nparticle, &
262 nparticle_kind, nparticle_local
263 INTEGER,
ALLOCATABLE,
DIMENSION(:, :) :: ascii
270 CALL force_env_get(force_env, subsys=subsys, para_env=para_env)
271 CALL cp_subsys_get(subsys, atomic_kinds=atomic_kinds, local_particles=local_particles, &
274 IF (
ASSOCIATED(local_particles%local_particle_set))
THEN
275 nparticle_kind = atomic_kinds%n_els
276 nparticle = particles%n_els
281 DO iparticle = 1, nparticle
282 DO iparticle_kind = 1, nparticle_kind
283 nparticle_local = local_particles%n_el(iparticle_kind)
284 DO iparticle_local = 1, nparticle_local
285 IF (iparticle == local_particles%list(iparticle_kind)%array(iparticle_local))
THEN
286 CALL local_particles%local_particle_set(iparticle_kind)% &
287 rng(iparticle_local)%stream%dump(rng_record=rng_record)
294 CALL para_env%sum(ascii)
301 END SUBROUTINE update_rng_particle
314 write_binary_restart_file)
318 LOGICAL,
INTENT(IN) :: skip_vel_section, &
319 write_binary_restart_file
321 CHARACTER(LEN=*),
PARAMETER :: routinen =
'update_subsys'
323 CHARACTER(LEN=default_string_length) :: coord_file_name, unit_str
324 INTEGER :: coord_file_format, handle, output_unit
325 INTEGER,
DIMENSION(:),
POINTER :: multiple_unit_cell
326 LOGICAL :: scale, use_ref_cell
327 REAL(kind=
dp) :: conv_factor
337 NULLIFY (work_section, core_particles, particles, shell_particles, &
338 subsys, cell, para_env, multipoles)
339 CALL timeset(routinen, handle)
340 CALL force_env_get(force_env, subsys=subsys, cell=cell, para_env=para_env)
342 CALL cp_subsys_get(subsys, particles=particles, molecules=molecules, &
343 shell_particles=shell_particles, core_particles=core_particles, &
344 multipoles=multipoles)
348 ALLOCATE (multiple_unit_cell(3))
349 multiple_unit_cell = 1
351 i_vals_ptr=multiple_unit_cell)
358 IF (.NOT. write_binary_restart_file)
THEN
360 i_val=coord_file_format)
363 c_val=coord_file_name)
365 IF (para_env%is_source())
THEN
366 CALL open_file(file_name=trim(adjustl(coord_file_name)), &
367 file_action=
"READWRITE", &
368 file_form=
"FORMATTED", &
369 file_position=
"REWIND", &
370 file_status=
"UNKNOWN", &
371 unit_number=output_unit)
372 CALL dump_coordinates_cp2k(particles, molecules, cell, conv_factor, &
373 output_unit=output_unit, &
374 core_or_shell=.false., &
375 scaled_coordinates=scale)
379 CALL section_coord_val_set(work_section, particles, molecules, conv_factor, scale, &
385 IF (.NOT. skip_vel_section)
THEN
386 IF (.NOT. write_binary_restart_file)
THEN
394 IF (.NOT. write_binary_restart_file)
THEN
395 IF (
ASSOCIATED(shell_particles))
THEN
400 CALL section_coord_val_set(work_section, shell_particles, molecules, &
401 conv_factor, scale, cell, shell=.true.)
402 IF (.NOT. skip_vel_section)
THEN
407 IF (
ASSOCIATED(core_particles))
THEN
412 CALL section_coord_val_set(work_section, core_particles, molecules, &
413 conv_factor, scale, cell, shell=.true.)
414 IF (.NOT. skip_vel_section)
THEN
424 CALL update_cell_section(cell, cell_section=work_section)
427 use_ref_cell = .false.
428 SELECT CASE (force_env%in_use)
430 CALL get_qs_env(force_env%qs_env, cell_ref=cell, use_ref_cell=use_ref_cell)
432 CALL eip_env_get(force_env%eip_env, cell_ref=cell, use_ref_cell=use_ref_cell)
434 CALL nnp_env_get(force_env%nnp_env, cell_ref=cell, use_ref_cell=use_ref_cell)
436 IF (use_ref_cell)
THEN
438 CALL update_cell_section(cell, cell_section=work_section)
442 IF (
ASSOCIATED(multipoles))
THEN
445 IF (
SIZE(work_section%values, 2) == 1)
EXIT
448 IF (
ASSOCIATED(multipoles%dipoles))
THEN
450 CALL update_dipoles_section(multipoles%dipoles, work_section)
452 IF (
ASSOCIATED(multipoles%quadrupoles))
THEN
454 CALL update_quadrupoles_section(multipoles%quadrupoles, work_section)
458 CALL timestop(handle)
468 SUBROUTINE update_cell_section(cell, cell_section)
473 CHARACTER(LEN=*),
PARAMETER :: routinen =
'update_cell_section'
476 INTEGER,
DIMENSION(:),
POINTER :: iwork
477 REAL(kind=
dp),
DIMENSION(:),
POINTER :: work
479 NULLIFY (work, iwork)
480 CALL timeset(routinen, handle)
484 work(1:3) = cell%hmat(1:3, 1)
489 work(1:3) = cell%hmat(1:3, 2)
494 work(1:3) = cell%hmat(1:3, 3)
506 CALL timestop(handle)
508 END SUBROUTINE update_cell_section
523 SUBROUTINE section_coord_val_set(coord_section, particles, molecules, conv_factor, &
529 REAL(kind=
dp) :: conv_factor
530 LOGICAL,
INTENT(IN) :: scale
532 LOGICAL,
INTENT(IN),
OPTIONAL :: shell
534 CHARACTER(LEN=*),
PARAMETER :: routinen =
'section_coord_val_set'
536 CHARACTER(LEN=2*default_string_length) :: line
537 CHARACTER(LEN=default_string_length) :: my_tag, name
538 INTEGER :: handle, ik, imol, irk, last_atom, nlist
539 LOGICAL :: ldum, molname_generated, my_shell
540 REAL(kind=
dp),
DIMENSION(3) :: r, s
544 TYPE(
val_type),
POINTER :: my_val, old_val
546 CALL timeset(routinen, handle)
548 NULLIFY (my_val, old_val, section, vals)
550 IF (
PRESENT(shell)) my_shell = shell
551 cpassert(
ASSOCIATED(coord_section))
552 cpassert(coord_section%ref_count > 0)
553 section => coord_section%section
556 CALL cp_abort(__location__, &
557 "section "//trim(section%name)//
" does not contain keyword "// &
561 IF (
SIZE(coord_section%values, 2) == 1)
EXIT
564 vals => coord_section%values(ik, 1)%list
566 IF (
ASSOCIATED(vals))
THEN
571 DO irk = 1, particles%n_els
574 s = particles%els(irk)%r
581 WRITE (unit=line, fmt=
"(T7,A,3ES25.16E3,1X,I0)") &
582 trim(adjustl(name)), s(1:3), particles%els(irk)%atom_index
588 new_pos => new_pos%rest
590 old_val => new_pos%first_el
592 new_pos%first_el => my_val
599 NULLIFY (new_pos%rest)
601 new_pos => new_pos%rest
606 IF (last_atom < irk)
THEN
608 molecule_now => molecules%els(imol)
610 CALL get_molecule_kind(molecule_now%molecule_kind, molname_generated=molname_generated, &
612 IF (molname_generated) my_tag =
""
616 s = particles%els(irk)%r
623 IF (len_trim(my_tag) > 0)
THEN
624 WRITE (unit=line, fmt=
"(T7,A,3ES25.16E3,1X,A,1X,I0)") &
625 trim(adjustl(name)), s(1:3), trim(my_tag), imol
627 WRITE (unit=line, fmt=
"(T7,A,3ES25.16E3)") &
628 trim(adjustl(name)), s(1:3)
636 new_pos => new_pos%rest
638 old_val => new_pos%first_el
640 new_pos%first_el => my_val
647 NULLIFY (new_pos%rest)
649 new_pos => new_pos%rest
656 coord_section%values(ik, 1)%list => vals
658 CALL timestop(handle)
660 END SUBROUTINE section_coord_val_set
670 SUBROUTINE update_dipoles_section(dipoles, dipoles_section)
672 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: dipoles
675 CHARACTER(LEN=*),
PARAMETER :: routinen =
'update_dipoles_section'
677 INTEGER :: handle, ik, irk, nlist, nloop
678 REAL(kind=
dp),
DIMENSION(:),
POINTER :: work
681 TYPE(
val_type),
POINTER :: my_val, old_val
683 CALL timeset(routinen, handle)
684 NULLIFY (my_val, old_val, section, vals)
685 cpassert(
ASSOCIATED(dipoles_section))
686 cpassert(dipoles_section%ref_count > 0)
687 section => dipoles_section%section
690 CALL cp_abort(__location__, &
691 "section "//trim(section%name)//
" does not contain keyword "// &
695 nloop =
SIZE(dipoles, 2)
697 IF (
SIZE(dipoles_section%values, 2) == 1)
EXIT
700 vals => dipoles_section%values(ik, 1)%list
702 IF (
ASSOCIATED(vals))
THEN
708 work = dipoles(1:3, irk)
715 new_pos => new_pos%rest
717 old_val => new_pos%first_el
719 new_pos%first_el => my_val
726 NULLIFY (new_pos%rest)
728 new_pos => new_pos%rest
733 dipoles_section%values(ik, 1)%list => vals
735 CALL timestop(handle)
737 END SUBROUTINE update_dipoles_section
747 SUBROUTINE update_quadrupoles_section(quadrupoles, quadrupoles_section)
749 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: quadrupoles
752 CHARACTER(LEN=*),
PARAMETER :: routinen =
'update_quadrupoles_section'
754 INTEGER :: handle, i, ik, ind, irk, j, nlist, nloop
755 REAL(kind=
dp),
DIMENSION(:),
POINTER :: work
758 TYPE(
val_type),
POINTER :: my_val, old_val
760 CALL timeset(routinen, handle)
761 NULLIFY (my_val, old_val, section, vals)
762 cpassert(
ASSOCIATED(quadrupoles_section))
763 cpassert(quadrupoles_section%ref_count > 0)
764 section => quadrupoles_section%section
767 CALL cp_abort(__location__, &
768 "section "//trim(section%name)//
" does not contain keyword "// &
772 nloop =
SIZE(quadrupoles, 2)
774 IF (
SIZE(quadrupoles_section%values, 2) == 1)
EXIT
777 vals => quadrupoles_section%values(ik, 1)%list
779 IF (
ASSOCIATED(vals))
THEN
789 work(ind) = quadrupoles(j, i, irk)
798 new_pos => new_pos%rest
800 old_val => new_pos%first_el
802 new_pos%first_el => my_val
809 NULLIFY (new_pos%rest)
811 new_pos => new_pos%rest
817 quadrupoles_section%values(ik, 1)%list => vals
819 CALL timestop(handle)
821 END SUBROUTINE update_quadrupoles_section
838 SUBROUTINE dump_coordinates_cp2k(particles, molecules, cell, conv_factor, &
839 output_unit, core_or_shell, &
845 REAL(kind=
dp),
INTENT(IN) :: conv_factor
846 INTEGER,
INTENT(IN) :: output_unit
847 LOGICAL,
INTENT(IN) :: core_or_shell, scaled_coordinates
849 CHARACTER(LEN=*),
PARAMETER :: routinen =
'dump_coordinates_cp2k'
851 CHARACTER(LEN=default_string_length) :: kind_name, tag
852 INTEGER :: handle, imolecule, iparticle, last_atom
853 LOGICAL :: ldum, molname_generated
854 REAL(kind=
dp),
DIMENSION(3) :: r, s
857 CALL timeset(routinen, handle)
859 cpassert(
ASSOCIATED(particles))
860 IF (.NOT. core_or_shell)
THEN
861 cpassert(
ASSOCIATED(molecules))
863 IF (scaled_coordinates)
THEN
864 cpassert(
ASSOCIATED(cell))
871 DO iparticle = 1, particles%n_els
872 CALL get_atomic_kind(particles%els(iparticle)%atomic_kind, name=kind_name)
873 IF (.NOT. core_or_shell)
THEN
874 IF (iparticle > last_atom)
THEN
875 imolecule = imolecule + 1
876 molecule => molecules%els(imolecule)
879 molname_generated=molname_generated, &
881 IF (molname_generated) tag =
""
886 IF (scaled_coordinates)
THEN
890 r(1:3) = particles%els(iparticle)%r(1:3)*conv_factor
892 IF (core_or_shell)
THEN
893 WRITE (unit=output_unit, fmt=
"(A,3ES25.16E3,1X,I0)") &
894 trim(adjustl(kind_name)), r(1:3), particles%els(iparticle)%atom_index
896 IF (len_trim(tag) > 0)
THEN
897 WRITE (unit=output_unit, fmt=
"(A,3ES25.16E3,1X,A,1X,I0)") &
898 trim(adjustl(kind_name)), r(1:3), trim(tag), imolecule
900 WRITE (unit=output_unit, fmt=
"(A,3ES25.16E3)") &
901 trim(adjustl(kind_name)), r(1:3)
906 CALL timestop(handle)
908 END SUBROUTINE dump_coordinates_cp2k
represent a simple array based list of the given type
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind(atomic_kind, fist_potential, element_symbol, name, mass, kind_number, natom, atom_list, rcov, rvdw, z, qeff, apol, cpol, mm_radius, shell, shell_active, damping)
Get attributes of an atomic kind.
Handles all functions related to the CELL.
subroutine, public real_to_scaled(s, r, cell)
Transform real to scaled cell coordinates. s=h_inv*r.
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
Utility routines to open and close files. Tracking of preconnections.
subroutine, public open_file(file_name, file_status, file_form, file_action, file_position, file_pad, unit_number, debug, skip_get_unit_number, file_access)
Opens the requested file using a free unit number.
subroutine, public close_file(unit_number, file_status, keep_preconnection)
Close an open file given by its logical unit number. Optionally, keep the file and unit preconnected.
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
real(kind=dp) function, public cp_unit_from_cp2k(value, unit_str, defaults, power)
converts from the internal cp2k units to the given unit
stores a lists of integer that are local to a processor. The idea is that these integers represent ob...
The environment for the empirical interatomic potential methods.
subroutine, public eip_env_get(eip_env, eip_model, eip_energy, eip_energy_var, eip_forces, coord_avg, coord_var, count, subsys, atomic_kind_set, particle_set, local_particles, molecule_kind_set, molecule_set, local_molecules, eip_input, force_env_input, cell, cell_ref, use_ref_cell, eip_kinetic_energy, eip_potential_energy, virial)
Returns various attributes of the eip environment.
Interface for the force calculations.
integer, parameter, public use_qmmm
subroutine, public multiple_fe_list(force_env_sections, root_section, i_force_eval, nforce_eval)
returns the order of the multiple force_env
integer, parameter, public use_eip_force
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
integer, parameter, public use_qs_force
integer, parameter, public use_nnp_force
Defines the basic variable types.
integer, parameter, public dp
integer, parameter, public default_string_length
Interface to the message passing library MPI.
Define the molecule kind structure types and the corresponding functionality.
subroutine, public get_molecule_kind(molecule_kind, atom_list, bond_list, bend_list, ub_list, impr_list, opbend_list, colv_list, fixd_list, g3x3_list, g4x6_list, vsite_list, torsion_list, shell_list, name, mass, charge, kind_number, natom, nbend, nbond, nub, nimpr, nopbend, nconstraint, nconstraint_fixd, nfixd, ncolv, ng3x3, ng4x6, nvsite, nfixd_restraint, ng3x3_restraint, ng4x6_restraint, nvsite_restraint, nrestraints, nmolecule, nsgf, nshell, ntorsion, molecule_list, nelectron, nelectron_alpha, nelectron_beta, bond_kind_set, bend_kind_set, ub_kind_set, impr_kind_set, opbend_kind_set, torsion_kind_set, molname_generated)
Get informations about a molecule kind.
represent a simple array based list of the given type
Define the data structure for the molecule information.
subroutine, public get_molecule(molecule, molecule_kind, lmi, lci, lg3x3, lg4x6, lcolv, first_atom, last_atom, first_shell, last_shell)
Get components from a molecule data set.
Multipole structure: for multipole (fixed and induced) in FF based MD.
Data types for neural network potentials.
subroutine, public nnp_env_get(nnp_env, nnp_forces, subsys, atomic_kind_set, particle_set, local_particles, molecule_kind_set, molecule_set, local_molecules, nnp_input, force_env_input, cell, cell_ref, use_ref_cell, nnp_potential_energy, virial)
Returns various attributes of the nnp environment.
Parallel (pseudo)random number generator (RNG) for multiple streams and substreams of random numbers.
integer, parameter, public rng_record_length
represent a simple array based list of the given type
logical function, public qmmm_ff_precond_only_qm(id1, id2, id3, id4, is_link)
This function handles the atom names and modifies the "_QM_" prefix, in order to find the parameters ...
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, ecoul_1c, rho0_s_rs, rho0_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, harris_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs)
Get the QUICKSTEP environment.
Utilities for string manipulations.
subroutine, public string_to_ascii(string, nascii)
Convert a string to sequence of integer numbers.
represent a list of objects
Type defining parameters related to the simulation cell.
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
represent a list of objects