All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
How to add plumed to an MD code

Plumed ships with scripts that can be used to add it to many of the standard MD packages.

Obviously though, if no patch is provided for the MD code you use then you will have to write one yourself. Plumed has been designed so that it can be added to an MD code either statically (i.e. as a collection of objects) or as a dynamic library. For technical reasons it is NOT possible to pack all the plumed objects as a static library and link the library (if you really want to know it: when static librearies are linked the linker typically discards all the objects which are not used; plumed is using an automatic registration process which is not compatible with this choice). Furthermore, unlike the previous version of plumed, the plumed source code is not compiled at the same time as your MD code is compiled. Instead plumed now has its own makefile and is compiled separately to the MD engines that rely on it. This makes plumed patches considerably simpler as now they only do two things:

Modifying your makefile

Once the plumed source code has been compiled a file src/Plumed.inc is generated. This file should be included in your MD code's makefile as it informs the code where to find all the plumed source. There are three possible ways to link PLUMED: static (linking the .o files directly), shared (linking the libplumed.so library) or runtime (which links only a wrapper to plumed, whereas the location of the remaining part on plumed, i.e. the libplumedKernel.so, can be specified at runtime). The relevant variables are

The libplumedKernel.so basically contains the same objects as the libplumed.so except for the wrapper.

The simplest approach is to take advantage of the three variants src/Plumed.inc.static , src/Plumed.inc.shared and src/Plumed.inc.runtime . Including one of those, it is sufficient to add to the linker line the macro $(PLUMED_LOAD).

The best way to patch your MD file is:

Now you can patch/unpatch your MD code, automatically linking PLUMED, as if it was an officially released interface (even though at this stage you do not have any functionality of plumed yet in your code - see below how to add it).

Note that after you have this done, PLUMED creates a file in $PLUMED_ROOT/patches/ that is generally named name-of-the-code.diff ( if your code is named name-of-the-code in the "save" phase). You can also manually create, in the same directory, a file called name-of-the-code.config so that you can use to perform five actions (totally optional). The first one can be used to do some heuristic check to verify that the patch is for the correct MD engine and that the MD source tree is in the correct state (this typically means that it hass been already configured). Then there are two functions that can do some modifications to some files that might vary considerable (as the Makefiles produced by the config procedure that change with the machine) and do some actions after the patching procedure is performed (a sanity check or a makedependence script to be rerun). The last two ones should revert the effect of the previous ones in such a way that it is possible to smoothly revert the patch. Try to keep them consistent so as the revert feature works correctly. This file is in bash and typically looks like this

function plumed_preliminary_test(){
# Just check of the code is the right one: no error means it is ok 
  grep -q  name-of-the-code Somefile 
}

function plumed_before_patch(){
# typically do something on the Makefile
}

function plumed_after_patch(){
# typically call some makedepends if necessary
}

function plumed_after_revert(){
# revert exacty the effect of plumed_before_patch
# typically undoes what was done on the Makefile
# additionally, one could add a makedepends
}

function plumed_before_revert(){
# revert exactly the effect of plumed_after_patch
# however, a makedepends script should always go in plumed_after_revert
}
Attention
Be careful with these scripts. You have access to a lot of flexibility (adding files, modifying existing ones arbitrarily etc) but you should always take care of the revertibility with "patch -r".

Language dependence

The next step is to add proper calls to plumed into the MD source code.

The interface between plumed and the MD codes is written is plain C so it is compatible with codes that are written in plain C, C++ and fortran. To use this interface one should first create a plumed object (plumed_create()), then sent to it messages (using plumed_cmd()) and deallocate it (plumed_finalize()) at the end of the simulation. All of these routines have C++ and FORTRAN equivalents. Notice that in C the plumed object is represented using a struct (plumed), in C++ it is represented using a class (PLMD::Plumed) and in FORTRAN using a string of 32 characters (CHARACTER(LEN=32)). Also notice that whenever passing strings from FORTRAN you should add to them "char(0)" as a string terminator to be sure that the string is properly interpreted (e.g. "pass-this-string"//char(0)).

As people might get confused by the idea of "plumed being an object", we also added the possibility of using a "generic global plumed instance", which can be accessed using routines with a g in the name (plumed_g_create(), plumed_g_cmd() and plumed_g_finalize()).

This interface is very general, will not change in future releases, and is fully contained in Plumed.h and Plumed.c files. For a reference to it, see Reference for interfacing MD codes with PLUMED

A practical example

We now describe how to pass data to plumed in from an MD code that is written in C/C++. First save all the files that you intend to modify into a .preplumed file (for example if you want to modify myfile.c save it first in its original version into myfile.c.preplumed): this will tell the patching procedure that this is a file that will enter the set of the files to be modified.

In C or C++ files containing calls to plumed you will have to include the Plumed.h file. In addition, you will also have to define a plumed obect that is visible in all these routines and you will probably like to define some sort of plumedswitch, which can be read in from the input to your code and used to tell the code whether or not this is to be a run with plumed. Finally, you might like to include something in input so that you specify the name of the plumed input file. How these things are best done will depend on your code and so we leave it to your discretion.

Plumed must perform three actions inside your MD code:

The only routine which can be used to send message to plumed is plumed_cmd (or, equivalently, Plumed::cmd in C++ and plumed_f_cmd in FORTRAN).

The various calls that can be used during initialization are as follows:

plumed plumedmain; plumedmain=plumed_create();                 // Create the plumed object

// Calls to pass data to plumed
plumed_cmd(plumedmain,"setRealPrecision",&real_precision);     // Pass a pointer to an integer containing the size of a real number (4 or 8)
plumed_cmd(plumedmain,"setMDEnergyUnits",&energyUnits);        // Pass a pointer to the conversion factor between the energy unit used in your code and kJ mol-1
plumed_cmd(plumedmain,"setMDLengthUnits",&lengthUnits);        // Pass a pointer to the conversion factor between the length unit used in your code and nm 
plumed_cmd(plumedmain,"setMDTimeUnits",&timeUnits);            // Pass a pointer to the conversion factor between the time unit used in your code and ps
plumed_cmd(plumedmain,"setPlumedDat",&plumedInput);            // Pass the name of the plumed input file from the md code to plumed
plumed_cmd(plumedmain,"setMPIComm",&MPI_COMM_WORLD);           // Pass a pointer to the MPI communicator to plumed
// notice that from fortran the command "setMPIFComm" should be used instead
plumed_cmd(plumedmain,"setNatoms",&natoms);                    // Pass a pointer to the number of atoms in the system to plumed
plumed_cmd(plumedmain,"setMDEngine","gromacs");                // Pass the name of your md engine to plumed (now it is just a label) 
plumed_cmd(plumedmain,"setLog",fplog);                         // Pass the file on which to write out the plumed log (if the file is already open)
plumed_cmd(plumedmain,"setLogFile",fplog);                     // Pass the file  on which to write out the plumed log (to be created)
plumed_cmd(plumedmain,"setTimestep",&delta_t);                 // Pass a pointer to the molecular dynamics timestep to plumed

// Calls to do the actual initialization (all the above commands must appear before this call)
plumed_cmd(plumedmain,"init",NULL);                            // Do all the initialization of plumed
plumed_cmd(plumedmain,"read",read);                            // Read the plumed input.  N.B. This is called during init and so this call is only required in special cases. 

Please note that if your code is in FORTRAN you should append a "char(0)" token to every string. Also, remember that FORTRAN is by default passing arguments by reference, so that the "&" symbols which are required in C are not necessary in FORTRAN.

The various calls that can be used pass data and calculate the forces due to the bias are as follows:

// Calls to pass data to plumed
plumed_cmd(plumedmain,"setStep",&step);                      // Pass a pointer to the current timestep to plumed
/ *** The way that you pass positions will depend on how they are stored in your code.  If the x, y and z position are all stored in a single array you may use:
plumed_cmd(plumedmain,"setPositions",&pos[0][0]);            // Pass a pointer to the first element in the atomic positions array to plumed  
                                                             // assuming they
                                                             // are stored in
                                                             // a
                                                             // x1,y1,z1,x2,y2,z2 ...
                                                             // kind of ordering
/ *** Othersize if you pass the three separate vectors of x, y and z positions using:
plumed_cmd(plumedmain,"setPositionX",&x[0]);                 // Pass a pointer to the first element in the array of x component of the atomic positions to plumed
plumed_cmd(plumedmain,"setPositionY",&y[0]);                 // Pass a pointer to the first element in the array of y component of the atomic positions to plumed
plumed_cmd(plumedmain,"setPositionZ",&z[0]);                 // Pass a pointer to the first element in the array of z component of the atomic positions to plumed
plumed_cmd(plumedmain,"setMasses",&mass[0]);                 // Pass a pointer to the first element in the masses array to plumed
plumed_cmd(plumedmain,"setCharges",&charge[0]);              // Pass a pointer to the first element in the charges array to plumed
plumed_cmd(plumedmain,"setBox",&box[0][0]);                  // Pass a pointer to the first element in the box share array to plumed
plumed_cmd(plumedmain,"setEnergy",&poteng);                  // Pass a pointer to the current value of the potential energy to plumed?
/ *** The way that you pass forces will depend on how they are stored in your code.  If the x, y and z force are all stored in a single array you may use:
plumed_cmd(plumedmain,"setForces",&f[0][0]);                 // Pass a pointer to the first element in the foces array to plumed
/ *** Othersize if you pass the three separate vectors of x, y and z forces using:
plumed_cmd(plumedmain,"setForcesX",&fx[0]);                  // Pass a pointer to the first element in the array of the x components of the atomic forces to plumed
plumed_cmd(plumedmain,"setForcesY",&fy[0]);                  // Pass a pointer to the first element in the array of the y components of the atomic forces to plumed
plumed_cmd(plumedmain,"setForcesZ",&fz[0]);                  // Pass a pointer to the first element in the array of the z components of the atomic forces to plumed
plumed_cmd(plumedmain,"setVirial",&force_vir[0][0]);         // Pass a pointer to the first element in the virial array to plumed

// Calls to do actual calculations
plumed_cmd(plumedmain,"calc",NULL);                          // Calculate and apply forces from the biases defined in the plumed input
// In some cases it can be faster to break up the above command into its individual setps:
plumed_cmd(plumedmain,"prepareCalc",NULL);                   // Prepare to do a calculation by requesting all the atomic positions from the MD code
plumed_cmd(plumedmain,"prepareDependencies",NULL);           // Work out what we are calculating during this MD step (this is the first step of prepareCalc)
plumed_cmd(plumedmain,"shareData",NULL);                     // Request all the atomic positions from the MD code (this is the second step of prepareCalc)
plumed_cmd(plumedmain,"performCalc",NULL);                   // Use the atomic positions collected during prepareCalc phase to calculate colvars and biases.

// Some extra calls that might come in handy
plumed_cmd(plumedmain,"createFullList",&n);                  // Create a list containing of all the atoms plumed is using to do calculations (return the number of atoms in n)
plumed_cmd(plumedmain,"getFullList",&list);                  // Return a list (in list) containing all the indices plumed is using to do calculations
plumed_cmd(plumedmain,"clearFullList",NULL);                 // Clear the list of all the atoms that plumed is using to do calculations
plumed_cmd(plumedmain,"clear",clear);                        // Clear and delete all the pointers inside plumed.

The plumed calls for the finalization tasks is as follows:

plumed_finalize(plumedmain);          // Call the plumed destructor

Dealing with parallelism

Plumed has functionality to deal with parallel MD codes. The particular form of the functionality used to do this (and the frequency with which you will have to call these routines) will depend on whether your code is parallelism using a domain decomposition or particle decomposition strategy. The calls for required for using this functionality are as follows:

plumed_cmd(plumedmain,"setAtomsNlocal",&nlocal);            // Pass a pointer to the number of atoms on this node
plumed_cmd(plumedmain,"setAtomsGatindex",gatindex);         // Pass an array containing the indices of all the atoms on this node (used for domain decomposition)
plumed_cmd(plumedmain,"setAtomsContiguous",&start);         // Number the atoms on this node from start to start+nlocal   (used for particle decomposition)

Inquiring for the plumed version

New functionalities might be added in the future to plumed. The description of all the possible "commands" sent to plumed is called its API (application programming interface). If you want to know which API is supported by plumed you can use this call:

plumed_cmd(plumedmain,"getApiVersion",&api);                 // Pass the api version that plumed is using

With the current version, this will set the api variable (an integer) to 1. As we add new features, this number will be increased.

the diffs

This is similar to plumed 1. All the files that you want to modify should be first copied to .preplumed files. Then use "plumed patch -s" to save the diff.

Notes

Compiler options to use shared libraries on many architectures: http://www.fortran-2000.com/ArnaudRecipes/sharedlib.html