Exemplo n.º 1
0
void ReferenceAtoms::readAtomsFromPDB( const PDB& pdb, const bool allowblocks  ){
  if( !allowblocks && pdb.getNumberOfAtomBlocks()!=1 ) error("found multi-atom-block pdb format but expecting only one block of atoms");  

  for(unsigned i=0;i<pdb.size();++i){
     indices.push_back( pdb.getAtomNumbers()[i] ); reference_atoms.push_back( pdb.getPositions()[i] );
     align.push_back( pdb.getOccupancy()[i] ); displace.push_back( pdb.getBeta()[i] );
  }
  atom_der_index.resize( reference_atoms.size() );
}
Exemplo n.º 2
0
PathMSDBase::PathMSDBase(const ActionOptions&ao):
PLUMED_COLVAR_INIT(ao),
neigh_size(-1),
neigh_stride(-1),
nframes(0)
{
  parse("LAMBDA",lambda);
  parse("NEIGH_SIZE",neigh_size);
  parse("NEIGH_STRIDE",neigh_stride);
  parse("REFERENCE",reference);

  // open the file
  FILE* fp=fopen(reference.c_str(),"r");
  std::vector<AtomNumber> aaa;
  if (fp!=NULL)
  {
    log<<"Opening reference file "<<reference.c_str()<<"\n";
    bool do_read=true;
    while (do_read){
         PDB mypdb; 
         RMSD mymsd(log); 
         do_read=mypdb.readFromFilepointer(fp,plumed.getAtoms().usingNaturalUnits(),0.1/atoms.getUnits().getLength());
         if(do_read){
            unsigned nat=0;
            nframes++;
            if(mypdb.getAtomNumbers().size()==0) error("number of atoms in a frame should be more than zero");
            if(nat==0) nat=mypdb.getAtomNumbers().size();
            if(nat!=mypdb.getAtomNumbers().size()) error("frames should have the same number of atoms");
            if(aaa.empty()) aaa=mypdb.getAtomNumbers();
            if(aaa!=mypdb.getAtomNumbers()) error("frames should contain same atoms in same order");
            log<<"Found PDB: "<<nframes<<" containing  "<<mypdb.getAtomNumbers().size()<<" atoms\n"; 
	    pdbv.push_back(mypdb); 
//            requestAtoms(mypdb.getAtomNumbers()); // is done in non base classes 
            derivs_s.resize(mypdb.getAtomNumbers().size());
            derivs_z.resize(mypdb.getAtomNumbers().size());
            mymsd.set(mypdb,"OPTIMAL");
            msdv.push_back(mymsd); // the vector that stores the frames
            //log<<mypdb; 
         }else{break ;}
    }
    fclose (fp);
    log<<"Found TOTAL "<<nframes<< " PDB in the file "<<reference.c_str()<<" \n"; 
    if(nframes==0) error("at least one frame expected");
  } 
  if(neigh_stride>0 || neigh_size>0){
           if(neigh_size>int(nframes)){
           	log.printf(" List size required ( %d ) is too large: resizing to the maximum number of frames required: %u  \n",neigh_size,nframes);
 		neigh_size=nframes;
           }
           log.printf("  Neighbor list enabled: \n");
           log.printf("                size   :  %d elements\n",neigh_size);
           log.printf("                stride :  %d timesteps \n",neigh_stride);
  }else{
           log.printf("  Neighbor list NOT enabled \n");
  }

}
Exemplo n.º 3
0
void ActionAtomistic::readAtomsFromPDB( const PDB& pdb ){
  Colvar*cc=dynamic_cast<Colvar*>(this);
  if(cc && cc->checkIsEnergy()) error("can't read energies from pdb files");

  for(unsigned j=0;j<indexes.size();j++){
      if( indexes[j].index()>pdb.size() ) error("there are not enough atoms in the input pdb file");
      if( pdb.getAtomNumbers()[j].index()!=indexes[j].index() ) error("there are atoms missing in the pdb file");  
      positions[j]=pdb.getPositions()[indexes[j].index()];
  }
  for(unsigned j=0;j<indexes.size();j++) charges[j]=pdb.getBeta()[indexes[j].index()];
  for(unsigned j=0;j<indexes.size();j++) masses[j]=pdb.getOccupancy()[indexes[j].index()];
}
bool DataCollectionObject::transferDataToPDB( PDB& mypdb ) {
  // Check if PDB contains argument names
  std::vector<std::string> pdb_args( mypdb.getArgumentNames() );
  // Now set the argument values
  std::map<std::string,double>::iterator it;
  for(unsigned i=0; i<pdb_args.size(); ++i) {
    it=args.find( pdb_args[i] );
    if( it==args.end() ) return false;
    mypdb.setArgumentValue( pdb_args[i], it->second );
  }
  // Now set the atomic positions
  std::vector<AtomNumber> pdb_pos( mypdb.getAtomNumbers() );
  if( pdb_pos.size()==positions.size() ) mypdb.setAtomPositions( positions );
  else if( pdb_pos.size()>0 ) plumed_merror("This feature is currently not ready");
  return true;
}
Exemplo n.º 5
0
void MultiDomainRMSD::read( const PDB& pdb ) {
  unsigned nblocks =  pdb.getNumberOfAtomBlocks();
  if( nblocks<2 ) error("multidomain RMSD only has one block of atoms");

  std::vector<Vector> positions; std::vector<double> align, displace;
  std::string num; blocks.resize( nblocks+1 ); blocks[0]=0;
  for(unsigned i=0; i<nblocks; ++i) blocks[i+1]=pdb.getAtomBlockEnds()[i];

  double tmp, lower=0.0, upper=std::numeric_limits<double>::max( );
  if( pdb.getArgumentValue("LOWER_CUTOFF",tmp) ) lower=tmp;
  if( pdb.getArgumentValue("UPPER_CUTOFF",tmp) ) upper=tmp;
  bool nopbc=pdb.hasFlag("NOPBC");

  domains.resize(0); weights.resize(0);
  for(unsigned i=1; i<=nblocks; ++i) {
    Tools::convert(i,num);
    if( ftype=="RMSD" ) {
      // parse("TYPE"+num, ftype );
      lower=0.0; upper=std::numeric_limits<double>::max( );
      if( pdb.getArgumentValue("LOWER_CUTOFF"+num,tmp) ) lower=tmp;
      if( pdb.getArgumentValue("UPPER_CUTOFF"+num,tmp) ) upper=tmp;
      nopbc=pdb.hasFlag("NOPBC");
    }
    domains.emplace_back( metricRegister().create<SingleDomainRMSD>( ftype ) );
    positions.resize( blocks[i] - blocks[i-1] );
    align.resize( blocks[i] - blocks[i-1] );
    displace.resize( blocks[i] - blocks[i-1] );
    unsigned n=0;
    for(unsigned j=blocks[i-1]; j<blocks[i]; ++j) {
      positions[n]=pdb.getPositions()[j];
      align[n]=pdb.getOccupancy()[j];
      displace[n]=pdb.getBeta()[j];
      n++;
    }
    domains[i-1]->setBoundsOnDistances( !nopbc, lower, upper );
    domains[i-1]->setReferenceAtoms( positions, align, displace );
    domains[i-1]->setupRMSDObject();

    double ww=0;
    if( !pdb.getArgumentValue("WEIGHT"+num,ww) ) weights.push_back( 1.0 );
    else weights.push_back( ww );
  }
  // And set the atom numbers for this object
  indices.resize(0); atom_der_index.resize(0);
  for(unsigned i=0; i<pdb.size(); ++i) { indices.push_back( pdb.getAtomNumbers()[i] ); atom_der_index.push_back(i); }
  // setAtomNumbers( pdb.getAtomNumbers() );
}
Exemplo n.º 6
0
void MultiDomainRMSD::read( const PDB& pdb ){
   unsigned nblocks =  pdb.getNumberOfAtomBlocks();
   if( nblocks<2 ) error("multidomain RMSD only has one block of atoms");
  
   std::vector<AtomNumber> atomnumbers;
   std::vector<Vector> positions; std::vector<double> align, displace;
   std::string num; blocks.resize( nblocks+1 ); blocks[0]=0;
   for(unsigned i=0;i<nblocks;++i) blocks[i+1]=pdb.getAtomBlockEnds()[i]; 

   double lower=0.0, upper=std::numeric_limits<double>::max( );
   parse("LOWER_CUTOFF",lower,true); 
   parse("UPPER_CUTOFF",upper,true);

   for(unsigned i=1;i<=nblocks;++i){
       Tools::convert(i,num);
       if( ftype=="RMSD" ){
          parse("TYPE"+num, ftype );
          parse("LOWER_CUTOFF"+num,lower,true); 
          parse("UPPER_CUTOFF"+num,upper,true); 
       }
       domains.push_back( metricRegister().create<SingleDomainRMSD>( ftype ) );
       positions.resize( blocks[i] - blocks[i-1] + 1 );
       align.resize( blocks[i] - blocks[i-1] + 1 );
       displace.resize( blocks[i] - blocks[i-1] + 1 );
       unsigned n=0;
       for(unsigned j=blocks[i-1];j<blocks[i];++j){
           positions[n]=pdb.getPositions()[j];
           align[n]=pdb.getOccupancy()[j];
           displace[n]=pdb.getBeta()[j];
           n++;
       }
       domains[i-1]->setBoundsOnDistances( true, lower, upper );  // Currently no option for nopbc
       domains[i-1]->setReferenceAtoms( positions, align, displace );
       domains[i-1]->setNumberOfAtoms( positions.size() );
       
       double ww=0; parse("WEIGHT"+num, ww, true );
       if( ww==0 ) weights.push_back( 1.0 );
       else weights.push_back( ww );
   }   
   // And set the atom numbers for this object
   setAtomNumbers( pdb.getAtomNumbers() );
}
Exemplo n.º 7
0
int main(int argc, char* argv[]) {

  // findiff step: eps
  double eps=1.e-6;

  // various setups: com needs to be reset
  bool remove_com=true;
  // normalize weights or not (to check if the rmsd is proportional to weights) 
  bool normalize_weights=true;
  // use msd instead of rmsd to check consistency
  bool squared=false;
  // enhance com to enphasize the issues with COM treatment
  bool enhance_com=false;

  // parse input instructions
  // default task, calculate the RMSD of one frame respect to a set of others 
  vector<int> task(1);task[0]=0;

  // this test wants to be only for OPTIMAL case
  string type; type.assign("OPTIMAL"); 

  // first parse the task: in this applications the tasks are simple integers that pass the action that is required
  for(int i = 1; i < argc; i++){ 
      task.push_back(atoi(argv[i]));
  }  
  if(std::find(task.begin(), task.end(), -1)!=task.end()){cout<<"squared=true (default false)"<<endl;squared=true;}
  if(std::find(task.begin(), task.end(), -2)!=task.end()){cout<<"normalize_weights=false (default true)"<<endl;normalize_weights=false;}
  if(std::find(task.begin(), task.end(), -3)!=task.end()){cout<<"remove_com=false (default true)"<<endl; remove_com=false;}
  if(std::find(task.begin(), task.end(), -4)!=task.end()){cout<<"OPTIMAL-FAST (default OPTIMAL)"<<endl; type.assign("OPTIMAL-FAST");}
  if(std::find(task.begin(), task.end(), -5)!=task.end()){cout<<"enhance_com=true (default false) include option -3 (no com removal) "<<endl;enhance_com=true ; }
  if(enhance_com)remove_com=false;


  cout<<"ARGUMENTS: \n";
  cout<<"OPTIONS that go on top of tasks:\n";
  cout<<" -1 : squared=true (default=false)\n";  
  cout<<" -2 : normalize_weights=false (default=true)\n";
  cout<<" -3 : remove_com=false (default=true) \n";
  cout<<" -4 : OPTIMAL-FAST (default=OPTIMAL) \n";
  cout<<" -5 : enhance_com=true (default false) automatically set option -3 (i.e. check the trouble with com when you do not remove it)\n";
  cout<<"TASKS (can choose more than one):\n";
  cout<<"  0 : normal rmsd/msd calculation  and derivative dumps (default: always done)\n";
  cout<<"  1 : findiff test for  d msd / d position  (inhomogenehous weights)\n";
  cout<<"  2 : findiff test for  d msd / d reference (inhomogenehous weights)\n";
  cout<<"  3 : findiff test for  d msd / d position  (homogenehous weights)\n";
  cout<<"  4 : findiff test for  d msd / d reference (homogenehous weights)\n";
  cout<<"  5 : findiff test for  d Rot / d position  (inhomogenehous weights)  (reference->position)\n";
  cout<<"  6 : findiff test for  d Rot / d reference  (inhomogenehous weights) (reference->position)\n";
  cout<<"  7 : consistency check for MSD proportionality (works with squared=true through option -1 )\n";
  cout<<"  8 : do some timings for all the above routines and for a growing number of atoms\n";
  cout<<"  9 : test the rotation order: print position.pdb reference.pdb aligned.pdb and check that it makes sense(should be reference aligned onto positions)\n";
  cout<<" 10 : findiff test for  d Rot / d position  ( position -> reference ) \n";
  cout<<" 11 : findiff test for  d Rot / d position  (homogenehous weights)  (reference->position)\n";
  cout<<" 12 : findiff test for  d Rot / d reference  (homogenehous weights) (reference->position)\n";
  cout<<" 13 : do timings only for the most common use (aligment +derivatives) and repeat for homogeneous weights\n";
  cout<<" 20 : a simple test calculating only derivative of positions which is needed to compare new version and old version (need to compile with -DOLDRMSD both plumed and the test)\n";

  PDB pdbref;

  // now create the object: does not do anything but set the typer to SIMPLE
  PLMD::RMSD* rmsd=new RMSD();

  // set the reference pdb 
  string reference; reference.assign("1GB1_mdl1_rototranslated.pdb");
  PDB pdb;
  if( !pdb.read(reference,false,0.1) ) 
        cout<<"missing input file 1GB1_mdl1_rototranslated.pdb "<<"\n";
  if (enhance_com){
	vector<Vector> v=pdb.getPositions();
	for(unsigned i=0;i<v.size();i++){v[i][0]+=10.;v[i][1]+=20.;v[i][2]+=30.;}
	pdb.setPositions(v);
  }

  
  cout<<"NOW CREATING THE RMSD OBJECT... with set() method";  
  // remember that "set" method parses the reference, align and displace from the PDB by calling the following methods 
  //      setReference(pdb.getPositions());  -> set the reference and put the weights as 1/n, remove the com according to such weight 
  //      setAlign(pdb.getOccupancy());    -> normalizes, set the alignment vector and remove the Com using such weights 
  //      setDisplace(pdb.getBeta());    -> normalizes and set the displacement vector 
  //      setType(mytype);
  rmsd->set(pdb,type,remove_com,normalize_weights);
  // store other vectors from further manipulation
  std::vector<Vector> ref ;  ref=pdb.getPositions() ;
  std::vector<double> align; align=pdb.getOccupancy(); // non-normalized !
  std::vector<double> displace;  displace=pdb.getBeta(); // non-normalized ! 

  cout<<"DONE!"<<endl;  

  // now take another conformation to compare with: the running frame
  PDB pdbrun; 
  // mimic gromacs: do it in nm
  if( !pdbrun.read("1GB1_mdl2.pdb",false,0.1) )
        cout<<"missing input file 1GB1_mdl2.pdb\n" ;
  std::vector<Vector> run ;  run=pdbrun.getPositions() ;
  std::vector<Vector> derivatives ; 
  if (enhance_com){
        for(unsigned i=0;i<run.size();i++){run[i][0]-=10.;run[i][1]-=20.;run[i][2]-=30.;}
  }




// Task 0: calculate the alignment and dump some data
  if(std::find(task.begin(), task.end(), 0)!=task.end()){

      cout<<"Task 0: calculates the alignment and retrieve some data"<<endl;
      double r=rmsd->calculate( run, derivatives, squared );
      cout<<"RMSD IS "<<r<<"\n";
      // now dump some more information
      ofstream myfile;
      myfile.open ("output_rmsd");
      myfile<<"RMSD ";
      myfile.setf( std::ios::fixed, std:: ios::floatfield );
      myfile.width(12);myfile.precision(6); myfile<<std::right<<r<<"\n";
      // the derivatives
      for(unsigned int i=0;i<run.size();++i){
       myfile<<"DDIST_DPOS "<<std::setw(5)<<i<<" X= "<<std::setw(18)<<std::right<<derivatives[i][0]<<" Y= "<<std::setw(18)<<std::right<<derivatives[i][1]<<" Z= "<<std::setw(18)<<std::right<<derivatives[i][2]<<"\n";
      }

      std::vector<Vector> DDistDRef; DDistDRef.resize(align.size());
      r=rmsd->calc_DDistDRef( run,derivatives,DDistDRef,squared);
      for(unsigned int i=0;i<run.size();++i){
       myfile<<"DDIST_DREF "<<std::setw(5)<<i<<" X= "<<std::setw(18)<<std::right<<DDistDRef[i][0]<<" Y= "<<std::setw(18)<<std::right<<DDistDRef[i][1]<<" Z= "<<std::setw(18)<<std::right<<DDistDRef[i][2]<<"\n";
      }
      myfile.close();
  }

  // Task 1: calculate findiff of running frame
  if(std::find(task.begin(), task.end(), 1)!=task.end()){
	cout<<"Task 1: calculates the finite difference for the running frame"<<endl;
  	double r_old=rmsd->calculate( run, derivatives, squared ); 
	std::vector<Vector> run_save=run;	
	for(unsigned int comp=0;comp<3;comp++){	
		for(unsigned int i=0;i<run.size();++i){
			//change the position		
			run[i][comp]+=eps;
  			double r=rmsd->calculate( run, derivatives, squared ); 
			cout<<"DDIST_DPOS COMPONENT "<<comp<<" "<<(r-r_old)/(run[i][comp]-run_save[i][comp])<<" "<<derivatives[i][comp]<<"\n";
			// restore the old position
			run=run_save;
		}
	}
  }

#ifndef OLDRMSD
  // Task 2: calculate findiff of reference frame
  if(std::find(task.begin(), task.end(), 2)!=task.end()){
	cout<<"Task 2: calculates the finite difference for the reference frame"<<endl;
  	double r_old=rmsd->calculate( run, derivatives, squared ); 
	std::vector<Vector> ref_save=ref;	
	std::vector<Vector> DDistDRef;	
        rmsd->clear();
        rmsd->set(align, displace, ref,type,remove_com,normalize_weights  );	
	for(unsigned int comp=0;comp<3;comp++){	
		for(unsigned int i=0;i<run.size();++i){
			//change the position		
			ref[i][comp]+=eps;
			// the asymmetry between reference and positions requires that 
			// in order to modify the reference one has to reset the rmsd object  
		        rmsd->clear();
			rmsd->set(align, displace, ref,type,remove_com,normalize_weights  );
  			double r=rmsd->calc_DDistDRef( run, derivatives, DDistDRef, squared); 
			cout<<"DDIST_DREF COMPONENT "<<comp<<" "<<(r-r_old)/(ref[i][comp]-ref_save[i][comp])<<" "<<DDistDRef[i][comp]<<"\n";
			// restore the old position
			ref=ref_save;
		}
	}
  }
 // Task 3 calculate findiff of running frame for alEqDis  version (align=displace)
  if(std::find(task.begin(), task.end(), 3)!=task.end()){
	std::vector<double> newalign(run.size(),1.); 
	std::vector<double> newdisplace(run.size(),1.); 
	rmsd->clear();
	rmsd->set(newalign, newdisplace, ref,type,remove_com,normalize_weights  );
	cout<<"Task 3: calculates the finite difference for the running frame (fast version)"<<endl;
  	double r_old=rmsd->calculate( run, derivatives, squared ); 
	std::vector<Vector> run_save=run;	
	for(unsigned int comp=0;comp<3;comp++){	
		for(unsigned int i=0;i<run.size();++i){
			//change the position		
			run[i][comp]+=eps;
  			double r=rmsd->calculate( run, derivatives, squared ); 
			cout<<"DDIST_DPOS_FAST COMPONENT "<<comp<<" "<<(r-r_old)/(run[i][comp]-run_save[i][comp])<<" "<<derivatives[i][comp]<<"\n";
			// restore the old position
			run=run_save;
		}
	}
  }

  // Task 4: calculate findiff of reference frame for alEqDis version (align=displace)
  if(std::find(task.begin(), task.end(), 4)!=task.end()){
	cout<<"Task 4: calculates the finite difference for the reference frame"<<endl;
	std::vector<double> newalign(run.size(),1.); 
	std::vector<double> newdisplace(run.size(),1.); 
	rmsd->clear();
	rmsd->set(newalign, newdisplace, ref,type,remove_com,normalize_weights);
  	double r_old=rmsd->calculate( run, derivatives, squared ); 
	std::vector<Vector> DDistDRef;	
	std::vector<Vector> ref_save=ref;	
	for(unsigned int comp=0;comp<3;comp++){	
		for(unsigned int i=0;i<run.size();++i){
			//change the position		
			ref[i][comp]+=eps;
			// this function below also reset the com of the reference (but not of the running frame)
			rmsd->clear();
			rmsd->set(newalign, newdisplace, ref,type,remove_com,normalize_weights  );
  			double r=rmsd->calc_DDistDRef( run, derivatives,DDistDRef, squared); 
			cout<<"DDIST_DREF_FAST COMPONENT "<<comp<<" "<<(r-r_old)/(ref[i][comp]-ref_save[i][comp])<<" "<<DDistDRef[i][comp]<<" "<<r<<" "<<r_old<<"\n";
			// restore the old position
			ref=ref_save;
		}
	}
  }
  
  // Task 5: calculate findiff of derivative of the rotation matrix respect to running frame
  if(std::find(task.begin(), task.end(), 5)!=task.end()){
	cout<<"Task 5: calculates the finite difference for derivative of the rotation matrix respect to the the running frame"<<endl;
        rmsd->clear();
        rmsd->set(align, displace, ref,type,remove_com,normalize_weights  );	

	Tensor Rotation,OldRotation;
	Matrix<std::vector<Vector> > DRotDPos(3,3);
        std::vector<Vector> DDistDRef;
	rmsd->calc_DDistDRef_Rot_DRotDPos( run, derivatives, DDistDRef, OldRotation, DRotDPos, squared ); 
	std::vector<Vector> run_save=run;	
	for(unsigned int a=0;a<3;a++){	
		for(unsigned int b=0;b<3;b++){	
			for(unsigned int comp=0;comp<3;comp++){	
				for(unsigned int i=0;i<run.size();++i){
					//change the position		
					run[i][comp]+=eps;
					rmsd->calc_DDistDRef_Rot_DRotDPos( run, derivatives ,DDistDRef, Rotation, DRotDPos, squared ); 
					cout<<"DROT_DPOS COMPONENT "<<comp<<" "<<(Rotation[a][b]-OldRotation[a][b])/(run[i][comp]-run_save[i][comp])<<" "<<DRotDPos[a][b][i][comp]<<"\n";
					// restore the old position
					run=run_save;
				}
			}
		}
	}
  }
  // Task 6: calculate findiff of derivative of the rotation matrix respect to reference frame 
  if(std::find(task.begin(), task.end(), 6)!=task.end()){
	cout<<"Task 6: calculates the finite difference for derivative of the rotation matrix respect to the the reference frame"<<endl;
        rmsd->clear();
        rmsd->set(align, displace, ref,type,remove_com,normalize_weights  );	
	Tensor Rotation,OldRotation;
	Matrix<std::vector<Vector> > DRotDPos(3,3),DRotDRef(3,3);
        std::vector<Vector> DDistDRef;
	rmsd->calc_DDistDRef_Rot_DRotDPos_DRotDRef( run, derivatives,  DDistDRef, OldRotation , DRotDPos , DRotDRef , squared); 
	std::vector<Vector> ref_save=ref;	
	for(unsigned int a=0;a<3;a++){	
		for(unsigned int b=0;b<3;b++){	
			for(unsigned int comp=0;comp<3;comp++){	
				for(unsigned int i=0;i<run.size();++i){
					//change the position		
					ref[i][comp]+=eps;
					// this function below also reset the com of the reference (but not of the running frame)
		                        rmsd->clear();
       			                rmsd->set(align, displace, ref,type,remove_com,normalize_weights  );
					rmsd->calc_DDistDRef_Rot_DRotDPos_DRotDRef( run, derivatives, DDistDRef, Rotation , DRotDPos, DRotDRef ,squared ); 
					cout<<"DROT_DREF COMPONENT "<<comp<<" "<<(Rotation[a][b]-OldRotation[a][b])/(ref[i][comp]-ref_save[i][comp])<<" "<<DRotDRef[a][b][i][comp]<<"\n";
					// restore the old position
					ref=ref_save;
				}
			}
		}
	}
  }
  // Task 7:  check weight consistency 

  if(std::find(task.begin(), task.end(), 7)!=task.end()){
	cout<<"Task 7: calculates the weight (displacement) consistency: all these should same result when weights are normalized in input by setReferenceAtoms otherwise they should be proportional when squared=true \n When squared=false, each factor of 2 in weights should produce a factor of sqrt(2) in the total value  "<<endl;
	rmsd->clear();
        rmsd->set(align, displace, ref,type,remove_com, false   );
  	double r=rmsd->calculate( run, derivatives, squared ); 
	cout<<"STANDARD WEIGHT "<<r<<"\n"; 

        std::vector<double> newdisplace=displace;for(std::vector<double>::iterator p=newdisplace.begin();p!=newdisplace.end();++p ){(*p)*=2.;} 
	rmsd->clear();
        rmsd->set(align, newdisplace, ref,type,remove_com, false   );
  	r=rmsd->calculate( run, derivatives,  squared ); 
	cout<<"DOUBLE WEIGHT "<<r<<"\n"; 

        newdisplace=displace;for(std::vector<double>::iterator p=newdisplace.begin();p!=newdisplace.end();++p ){(*p)*=4.;} 
	rmsd->clear();
        rmsd->set(align,newdisplace, ref,type,remove_com, false );
  	r=rmsd->calculate( run, derivatives, squared ); 
	cout<<"FOUR WEIGHT "<<r<<"\n"; 
  }

  // Task 8: do some timings to get a flavor
  if(std::find(task.begin(), task.end(), 8)!=task.end()){
      cout<<"Task 8: makes some timings for increasing atoms and different routines "<<endl;
      rmsd->clear();
      rmsd->set(align, displace, ref,type,remove_com,normalize_weights  );	
      vector<Vector> r_run,r_ref;
      vector<double> r_al,r_disp;
      for (unsigned int i=0;i<10;i++){r_run.push_back(run[i]);r_ref.push_back(ref[i]);r_al.push_back(align[i]);r_disp.push_back(displace[i]);}

      for(unsigned int i=0;i<10;i++){
      	cout<<"NUMBER OF ATOMS : "<<r_run.size()<<endl;
      	unsigned ntest; ntest=100;
      	// test the fast routine
	rmsd->clear();
        rmsd->set(r_al,r_disp, r_ref, type,remove_com, normalize_weights );
 
      	Stopwatch sw;
      	sw.start();	
        for(unsigned int j=0;j<ntest;j++)rmsd->calculate( r_run, derivatives, squared );
      	sw.stop();	
      	cout<<"SIMPLE ROUTINE \n"<<sw<<endl;

        std::vector<Vector> DDistDRef;
      	Stopwatch sw2;
      	sw2.start();	
        for(unsigned int j=0;j<ntest;j++)rmsd->calc_DDistDRef( r_run, derivatives,DDistDRef, squared); 
      	sw2.stop();	
      	cout<<"WITH REFERENCE FRAME: \n"<<sw2<<endl;

      	Tensor Rotation;
      	Matrix<std::vector<Vector> > DRotDPos(3,3);	
      	Stopwatch sw3;
      	sw3.start();	
        for(unsigned int j=0;j<ntest;j++)rmsd->calc_DDistDRef_Rot_DRotDPos( r_run, derivatives,DDistDRef, Rotation , DRotDPos, squared); 
      	sw3.stop();	
      	cout<<"WITH ROTATION MATRIX DERIVATIVE: \n"<<sw3<<endl;

        Matrix<std::vector<Vector> > DRotDRef(3,3);
      	Stopwatch sw4;
      	sw4.start();	
        for(unsigned int j=0;j<ntest;j++)rmsd->calc_DDistDRef_Rot_DRotDPos_DRotDRef( r_run, derivatives ,DDistDRef, Rotation , DRotDPos, DRotDRef, squared); 
      	sw4.stop();	
      	cout<<"WITH ROTATION MATRIX DERIVATIVE OF REEFERENCE: \n"<<sw4<<endl;
      	// duplicate the atoms
      	unsigned s=r_run.size();
      	for (unsigned int i=0;i<s;i++){r_run.push_back(r_run[i]);r_ref.push_back(r_ref[i]);r_al.push_back(r_al[i]);r_disp.push_back(r_disp[i]);}
      
      }
  } 
  // Task 9: check the rotation
  if(std::find(task.begin(), task.end(), 9)!=task.end()){
      cout<<"Task 9: dump some pdbs so to check if all makes sense when using inverse transform. In particular positions_aligned.pdb should overlap with reference_centered.pdb "<<endl;
	rmsd->clear();
        rmsd->set(align, displace, ref,type,remove_com,normalize_weights  );	
	// dump the reference
	ofstream myfile;
	myfile.open ("reference.pdb");		
	std::vector<AtomNumber> at=pdb.getAtomNumbers();
	std::vector<Vector>   pos=pdb.getPositions();
	unsigned k=0;
	for(std::vector<AtomNumber>::iterator i=at.begin(); i!=at.end(); i++){
		myfile<<"ATOM";
                myfile.width(7);myfile<<std::right<<(*i).serial()<<" "; 
                myfile.width(4);myfile<<std::left<<pdb.getAtomName(*i); 
                myfile.width(4);myfile<<std::right<<pdb.getResidueName(*i)<<" A"; 
                myfile.width(4);myfile<<std::right<<pdb.getResidueNumber(*i)<<"    "; 
		myfile.setf( std::ios::fixed, std:: ios::floatfield );
                myfile.width(8);myfile.precision(3); myfile<<std::right<<pos[k][0]*10; 
                myfile.width(8);myfile.precision(3); myfile<<std::right<<pos[k][1]*10; 
                myfile.width(8);myfile.precision(3); myfile<<std::right<<pos[k][2]*10<<"  1.00  1.00\n"; 
		k++;	
	}	
	myfile.close();			
	// dump the position
	myfile.open ("positions.pdb");		
	at=pdbrun.getAtomNumbers();
	std::vector<Vector>   runpos=pdbrun.getPositions();
	k=0;
	for(std::vector<AtomNumber>::iterator i=at.begin(); i!=at.end(); i++){
		myfile<<"ATOM";
                myfile.width(7);myfile<<std::right<<(*i).serial()<<" "; 
                myfile.width(4);myfile<<std::left<<pdbrun.getAtomName(*i); 
                myfile.width(4);myfile<<std::right<<pdbrun.getResidueName(*i)<<" A"; 
                myfile.width(4);myfile<<std::right<<pdbrun.getResidueNumber(*i)<<"    "; 
		myfile.setf( std::ios::fixed, std:: ios::floatfield );
                myfile.width(8);myfile.precision(3); myfile<<std::right<<runpos[k][0]*10; 
                myfile.width(8);myfile.precision(3); myfile<<std::right<<runpos[k][1]*10; 
                myfile.width(8);myfile.precision(3); myfile<<std::right<<runpos[k][2]*10<<"  1.00  1.00\n"; 
		k++;	
	}	
	myfile.close();			
	// now do the alignment
	Tensor Rotation;
	Matrix<std::vector<Vector> > DRotDPos(3,3);
        std::vector<Vector> DDistDRef;
	std::vector<Vector> alignedpos;
	std::vector<Vector> centeredpos;
	std::vector<Vector> centeredref;
	std::vector<Vector> ddistdpos;
	rmsd->calc_PCAelements( run, derivatives, Rotation ,  DRotDPos , alignedpos ,centeredpos, centeredref ,squared); 
	myfile.open ("positions_aligned.pdb");		
	k=0;
	for(std::vector<AtomNumber>::iterator i=at.begin(); i!=at.end(); i++){
		myfile<<"ATOM";
                myfile.width(7);myfile<<std::right<<(*i).serial()<<" "; 
                myfile.width(4);myfile<<std::left<<pdbrun.getAtomName(*i); 
                myfile.width(4);myfile<<std::right<<pdbrun.getResidueName(*i)<<" A"; 
                myfile.width(4);myfile<<std::right<<pdbrun.getResidueNumber(*i)<<"    "; 
		myfile.setf( std::ios::fixed, std:: ios::floatfield );
                myfile.width(8);myfile.precision(3); myfile<<std::right<<alignedpos[k][0]*10; 
                myfile.width(8);myfile.precision(3); myfile<<std::right<<alignedpos[k][1]*10; 
                myfile.width(8);myfile.precision(3); myfile<<std::right<<alignedpos[k][2]*10<<"  1.00  1.00\n"; 
		k++;	
	}	
	myfile.close();			
	// dump the aligned	
	myfile.open ("reference_centered.pdb");		
	k=0;
	for(std::vector<AtomNumber>::iterator i=at.begin(); i!=at.end(); i++){
		myfile<<"ATOM";
                myfile.width(7);myfile<<std::right<<(*i).serial()<<" "; 
                myfile.width(4);myfile<<std::left<<pdbrun.getAtomName(*i); 
                myfile.width(4);myfile<<std::right<<pdbrun.getResidueName(*i)<<" A"; 
                myfile.width(4);myfile<<std::right<<pdbrun.getResidueNumber(*i)<<"    "; 
		myfile.setf( std::ios::fixed, std:: ios::floatfield );
                myfile.width(8);myfile.precision(3); myfile<<std::right<<centeredref[k][0]*10; 
                myfile.width(8);myfile.precision(3); myfile<<std::right<<centeredref[k][1]*10; 
                myfile.width(8);myfile.precision(3); myfile<<std::right<<centeredref[k][2]*10<<"  1.00  1.00\n"; 
		k++;	
	}	
	myfile.close();			
  }
  // Task 10: derivative of the rotation matrix (in case of reverse transition) 
  if(std::find(task.begin(), task.end(), 10)!=task.end()){
	cout<<"Task 10: calculates the finite difference for derivative of the rotation matrix respect to the the running frame"<<endl;
	rmsd->clear();
        rmsd->set(align, displace, ref,type,remove_com,normalize_weights  );	
	Tensor Rotation,OldRotation;
	Matrix<std::vector<Vector> > DRotDPos(3,3);
        std::vector<Vector> DDistDRef;
        std::vector<Vector> alignedpos;
        std::vector<Vector> centeredpos;
        std::vector<Vector> centeredref;
	std::vector<Vector> ddistdpos;
	rmsd->calc_PCAelements( run, derivatives, OldRotation , DRotDPos , alignedpos ,centeredpos, centeredref, squared ); 
	std::vector<Vector> run_save=run;	
	for(unsigned int a=0;a<3;a++){	
		for(unsigned int b=0;b<3;b++){	
			for(unsigned int comp=0;comp<3;comp++){	
				for(unsigned int i=0;i<run.size();++i){
					//change the position		
					run[i][comp]+=eps;
					rmsd->calc_PCAelements( run, derivatives, Rotation ,DRotDPos , alignedpos ,centeredpos, centeredref , squared); 
					cout<<"DROT_DPOS_INVERSE_TRANSFORM COMPONENT "<<comp<<" "<<(Rotation[a][b]-OldRotation[a][b])/(run[i][comp]-run_save[i][comp])<<" "<<DRotDPos[a][b][i][comp]<<"\n";
					// restore the old position
					run=run_save;
				}
			}
		}
	}
  }
  // Task 11: calculate findiff of derivative of the rotation matrix respect to running frame (homogenehous weights)
  if(std::find(task.begin(), task.end(), 11)!=task.end()){
	cout<<"Task 11: calculates the finite difference for derivative of the rotation matrix respect to the the running frame (homogeneous weights)"<<endl;
	std::vector<double> newalign(run.size(),1.); 
	std::vector<double> newdisplace(run.size(),1.); 
	rmsd->clear();
        rmsd->set(newalign, newdisplace, ref,type,remove_com,normalize_weights  );	
	Tensor Rotation,OldRotation;
	Matrix<std::vector<Vector> > DRotDPos(3,3);
        std::vector<Vector> DDistDRef;
	std::vector<Vector> alignedpos;
	std::vector<Vector> centeredpos;
	std::vector<Vector> centeredref;
	rmsd->calc_DDistDRef_Rot_DRotDPos( run,derivatives, DDistDRef, OldRotation , DRotDPos, squared  ); 
	//rmsd->calc_PCAelements( run, derivatives, OldRotation ,  DRotDPos , alignedpos ,centeredpos, centeredref ,squared); 
	std::vector<Vector> run_save=run;	
	for(unsigned int a=0;a<3;a++){	
		for(unsigned int b=0;b<3;b++){	
			for(unsigned int comp=0;comp<3;comp++){	
				for(unsigned int i=0;i<run.size();++i){
					//change the position		
					run[i][comp]+=eps;
					rmsd->calc_DDistDRef_Rot_DRotDPos( run, derivatives , DDistDRef, Rotation , DRotDPos, squared  ); 
					//rmsd->calc_PCAelements( run, derivatives, Rotation ,  DRotDPos , alignedpos ,centeredpos, centeredref ,squared); 
					cout<<"DROT_DPOS COMPONENT "<<comp<<" "<<(Rotation[a][b]-OldRotation[a][b])/(run[i][comp]-run_save[i][comp])<<" "<<DRotDPos[a][b][i][comp]<<"\n";
					// restore the old position
					run=run_save;
				}
			}
		}
	}
  }
  // Task 12: calculate findiff of derivative of the rotation matrix respect to reference frame 
  if(std::find(task.begin(), task.end(), 12)!=task.end()){
	cout<<"Task 12: calculates the finite difference for derivative of the rotation matrix respect to the the reference frame (homogeneous weights)"<<endl;
	std::vector<double> newalign(run.size(),1.); 
	std::vector<double> newdisplace(run.size(),1.); 
	rmsd->clear();
        rmsd->set(newalign, newdisplace, ref,type,remove_com,normalize_weights  );	
	Tensor Rotation,OldRotation;
	Matrix<std::vector<Vector> > DRotDPos(3,3),DRotDRef(3,3);
        std::vector<Vector> DDistDRef;
	rmsd->calc_DDistDRef_Rot_DRotDPos_DRotDRef( run, derivatives, DDistDRef, OldRotation , DRotDPos , DRotDRef ,squared); 
	std::vector<Vector> ref_save=ref;	
	for(unsigned int a=0;a<3;a++){	
		for(unsigned int b=0;b<3;b++){	
			for(unsigned int comp=0;comp<3;comp++){	
				for(unsigned int i=0;i<run.size();++i){
					//change the position		
					ref[i][comp]+=eps;
					// this function below also reset the com of the reference (but not of the running frame)
					rmsd->clear();
				        rmsd->set(newalign, newdisplace, ref,type,remove_com,normalize_weights  );	
					rmsd->calc_DDistDRef_Rot_DRotDPos_DRotDRef( run, derivatives, DDistDRef, Rotation , DRotDPos, DRotDRef ,squared ); 
					cout<<"DROT_DREF COMPONENT "<<comp<<" "<<(Rotation[a][b]-OldRotation[a][b])/(ref[i][comp]-ref_save[i][comp])<<" "<<DRotDRef[a][b][i][comp]<<"\n";
					// restore the old position
					ref=ref_save;
				}
			}
		}
	}
  }
   // Task 13: do some timings to get a flavor (only on alignment+derivatives)
  if(std::find(task.begin(), task.end(), 13)!=task.end()){
	cout<<"Task 13: timings for the most common use (only on alignment+derivatives)  "<<endl;
	vector<Vector> r_run,r_ref;
	vector<double> r_al,r_disp;
	for (unsigned int i=0;i<10;i++){r_run.push_back(run[i]);r_ref.push_back(ref[i]);r_al.push_back(align[i]);r_disp.push_back(displace[i]);}
	for(unsigned int i=0;i<10;i++){
		cout<<"NUMBER OF ATOMS : "<<r_run.size()<<endl;
		unsigned ntest; ntest=100;
		// test the fast routine
		rmsd->clear();
	        rmsd->set(r_al, r_disp, r_ref, type,remove_com,normalize_weights  );	
		Stopwatch sw;
		sw.start();	
	        for(unsigned int j=0;j<ntest;j++)rmsd->calculate( r_run, derivatives, squared );
		sw.stop();	
                cout<<"TIME \n"<<sw<<endl;
		// duplicate the atoms
		unsigned s=r_run.size();
		for (unsigned int i=0;i<s;i++){r_run.push_back(r_run[i]);r_ref.push_back(r_ref[i]);r_al.push_back(r_al[i]);r_disp.push_back(r_disp[i]);}
	
	}
	cout <<"NOW HOMOGENEOUS WEIGHTS\n";
	std::vector<double> newalign(10,1.); 
	std::vector<double> newdisplace(10,1.); 
	r_run.resize(0);
	r_ref.resize(0);
	r_al.resize(0);
	r_disp.resize(0);
	for (unsigned int i=0;i<10;i++){r_run.push_back(run[i]);r_ref.push_back(ref[i]);r_al.push_back(newalign[i]);r_disp.push_back(newdisplace[i]);}
        rmsd->clear();
        rmsd->set(r_al,r_disp , r_ref ,type,remove_com,normalize_weights  );
	for(unsigned int i=0;i<10;i++){
		cout<<"NUMBER OF ATOMS : "<<r_run.size()<<endl;
		unsigned ntest; ntest=100;
		// test the fast routine
                rmsd->clear();
                rmsd->set(r_al, r_disp, r_ref, type,remove_com,normalize_weights  );
		derivatives.resize(r_al.size());
		Stopwatch sw;
		sw.start();	
	        for(unsigned int j=0;j<ntest;j++)rmsd->calculate( r_run, derivatives, squared );
		sw.stop();	
                cout<<"TIME \n"<<sw<<endl;
		// duplicate the atoms
		unsigned s=r_run.size();
		for (unsigned int i=0;i<s;i++){r_run.push_back(r_run[i]);r_ref.push_back(r_ref[i]);r_al.push_back(r_al[i]);r_disp.push_back(r_disp[i]);}
	
	}
	
	

  } 


#endif
  // Task 20: do some timings to get a flavor
  if(std::find(task.begin(), task.end(), 20)!=task.end()){
      cout<<"Task 8: makes some timings for increasing atoms to compare new and old rmsd (need to be recompiled with -DOLDRMSD) "<<endl;
      vector<Vector> r_run,r_ref;
      vector<double> r_al,r_disp;
      for (unsigned int i=0;i<10;i++){r_run.push_back(run[i]);r_ref.push_back(ref[i]);r_al.push_back(align[i]);r_disp.push_back(displace[i]);}

      for(unsigned int i=0;i<10;i++){
      	cout<<"NUMBER OF ATOMS : "<<r_run.size()<<endl;
      	unsigned ntest; ntest=100;
      	// test the fast routine
	rmsd->clear();
        rmsd->set(r_al,r_disp, r_ref, type,remove_com, normalize_weights );
 
      	Stopwatch sw;
      	sw.start();	
        for(unsigned int j=0;j<ntest;j++)rmsd->calculate( r_run, derivatives, squared );
      	sw.stop();	
      	cout<<"SIMPLE ROUTINE \n"<<sw<<endl;

      	unsigned s=r_run.size();
      	for (unsigned int i=0;i<s;i++){r_run.push_back(r_run[i]);r_ref.push_back(r_ref[i]);r_al.push_back(r_al[i]);r_disp.push_back(r_disp[i]);}
      
      }
  } 


  return 0;
}
Exemplo n.º 8
0
int Driver<real>::main(FILE* in,FILE*out,Communicator& pc){

  Units units;
  PDB pdb;

// Parse everything
  bool printhelpdebug; parseFlag("--help-debug",printhelpdebug);
  if( printhelpdebug ){
      fprintf(out,"%s",
         "Additional options for debug (only to be used in regtest):\n"
         "  [--debug-float]         : turns on the single precision version (to check float interface)\n"
         "  [--debug-dd]            : use a fake domain decomposition\n"
         "  [--debug-pd]            : use a fake particle decomposition\n"
      );
      return 0;
  }
  // Are we reading trajectory data
  bool noatoms; parseFlag("--noatoms",noatoms);

  std::string fakein; 
  bool debugfloat=parse("--debug-float",fakein);
  if(debugfloat && sizeof(real)!=sizeof(float)){
      CLTool* cl=cltoolRegister().create(CLToolOptions("driver-float"));    //new Driver<float>(*this);
      cl->setInputData(this->getInputData());
      int ret=cl->main(in,out,pc);
      delete cl;
      return ret;
  }

  bool debug_pd=parse("--debug-pd",fakein);
  bool debug_dd=parse("--debug-dd",fakein);
  if( debug_pd || debug_dd ){
    if(noatoms) error("cannot debug without atoms");
  }

// set up for multi replica driver:
  int multi=0;
  parse("--multi",multi);
  Communicator intracomm;
  Communicator intercomm;
  if(multi){
    int ntot=pc.Get_size();
    int nintra=ntot/multi;
    if(multi*nintra!=ntot) error("invalid number of processes for multi environment");
    pc.Split(pc.Get_rank()/nintra,pc.Get_rank(),intracomm);
    pc.Split(pc.Get_rank()%nintra,pc.Get_rank(),intercomm);
  } else {
    intracomm.Set_comm(pc.Get_comm());
  }

// set up for debug replica exchange:
  bool debug_grex=parse("--debug-grex",fakein);
  int  grex_stride=0;
  FILE*grex_log=NULL;
  if(debug_grex){
    if(noatoms) error("must have atoms to debug_grex");
    if(multi<2)  error("--debug_grex needs --multi with at least two replicas");
    Tools::convert(fakein,grex_stride);
    string n; Tools::convert(intercomm.Get_rank(),n);
    string file;
    parse("--debug-grex-log",file);
    if(file.length()>0){
      file+="."+n;
      grex_log=fopen(file.c_str(),"w");
    }
  }

// Read the plumed input file name  
  string plumedFile; parse("--plumed",plumedFile);
// the timestep
  double t; parse("--timestep",t);
  real timestep=real(t);
// the stride
  unsigned stride; parse("--trajectory-stride",stride);
// are we writing forces
  string dumpforces(""), dumpforcesFmt("%f");; 
  if(!noatoms) parse("--dump-forces",dumpforces);
  if(dumpforces!="") parse("--dump-forces-fmt",dumpforcesFmt);

  string trajectory_fmt;

// Read in an xyz file
  string trajectoryFile(""), pdbfile("");
  bool pbc_cli_given=false; vector<double> pbc_cli_box(9,0.0);
  if(!noatoms){
     std::string traj_xyz; parse("--ixyz",traj_xyz);
     std::string traj_gro; parse("--igro",traj_gro);
     if(traj_xyz.length()>0 && traj_gro.length()>0){
       fprintf(stderr,"ERROR: cannot provide more than one trajectory file\n");
       if(grex_log)fclose(grex_log);
       return 1;
     }
     if(traj_xyz.length()>0 && trajectoryFile.length()==0){
       trajectoryFile=traj_xyz;
       trajectory_fmt="xyz";
     }
     if(traj_gro.length()>0 && trajectoryFile.length()==0){
       trajectoryFile=traj_gro;
       trajectory_fmt="gro";
     }
     if(trajectoryFile.length()==0){
       fprintf(stderr,"ERROR: missing trajectory data\n"); 
       if(grex_log)fclose(grex_log);
       return 1;
     }
     string lengthUnits(""); parse("--length-units",lengthUnits);
     if(lengthUnits.length()>0) units.setLength(lengthUnits);
  
     parse("--pdb",pdbfile);
     if(pdbfile.length()>0){
       bool check=pdb.read(pdbfile,false,1.0);
       if(!check) error("error reading pdb file");
     }

     string pbc_cli_list; parse("--box",pbc_cli_list);
     if(pbc_cli_list.length()>0) {
       pbc_cli_given=true;
       vector<string> words=Tools::getWords(pbc_cli_list,",");
       if(words.size()==3){
         for(int i=0;i<3;i++) sscanf(words[i].c_str(),"%100lf",&(pbc_cli_box[4*i]));
       } else if(words.size()==9) {
         for(int i=0;i<9;i++) sscanf(words[i].c_str(),"%100lf",&(pbc_cli_box[i]));
       } else {
         string msg="ERROR: cannot parse command-line box "+pbc_cli_list;
         fprintf(stderr,"%s\n",msg.c_str());
         return 1;
       }

     }
  }

  if( debug_dd && debug_pd ) error("cannot use debug-dd and debug-pd at the same time");
  if(debug_pd || debug_dd){
    if( !Communicator::initialized() ) error("needs mpi for debug-pd");
  }

  Plumed p;
  int rr=sizeof(real);
  p.cmd("setRealPrecision",&rr);
  int checknatoms=-1;
  int step=0;
  if(Communicator::initialized()){
    if(multi){
      if(intracomm.Get_rank()==0) p.cmd("GREX setMPIIntercomm",&intercomm.Get_comm());
      p.cmd("GREX setMPIIntracomm",&intracomm.Get_comm());
      p.cmd("GREX init");
    } 
    p.cmd("setMPIComm",&intracomm.Get_comm());
  } 
  p.cmd("setMDLengthUnits",&units.getLength());
  p.cmd("setMDEngine","driver");
  p.cmd("setTimestep",&timestep);
  p.cmd("setPlumedDat",plumedFile.c_str());
  p.cmd("setLog",out);

  if(multi){
    string n;
    Tools::convert(intercomm.Get_rank(),n);
    trajectoryFile+="."+n;
  }

  FILE* fp=NULL; FILE* fp_forces=NULL;
  if(!noatoms){
     if (trajectoryFile=="-") 
       fp=in;
     else {
       fp=fopen(trajectoryFile.c_str(),"r");
       if(!fp){
         string msg="ERROR: Error opening XYZ file "+trajectoryFile;
         fprintf(stderr,"%s\n",msg.c_str());
         return 1;
       }
     }

     if(dumpforces.length()>0){
       if(Communicator::initialized() && pc.Get_size()>1){
         string n;
         Tools::convert(pc.Get_rank(),n);
         dumpforces+="."+n;
       }
       fp_forces=fopen(dumpforces.c_str(),"w");
     }
  }

  std::string line;
  std::vector<real> coordinates;
  std::vector<real> forces;
  std::vector<real> masses;
  std::vector<real> charges;
  std::vector<real> cell;
  std::vector<real> virial;

// variables to test particle decomposition
  int pd_nlocal;
  int pd_start;
// variables to test random decomposition (=domain decomposition)
  std::vector<int>  dd_gatindex;
  std::vector<int>  dd_g2l;
  std::vector<real> dd_masses;
  std::vector<real> dd_charges;
  std::vector<real> dd_forces;
  std::vector<real> dd_coordinates;
  int dd_nlocal;
// random stream to choose decompositions
  Random rnd;

  while(true){
    if(!noatoms){
       if(!Tools::getline(fp,line)) break;
    }

    int natoms;
    bool first_step=false;
    if(!noatoms){
      if(trajectory_fmt=="gro") if(!Tools::getline(fp,line)) error("premature end of trajectory file");
      sscanf(line.c_str(),"%100d",&natoms);
    }
    if(checknatoms<0 && !noatoms){
      pd_nlocal=natoms;
      pd_start=0;
      first_step=true;
      masses.assign(natoms,real(1.0));
      charges.assign(natoms,real(0.0));
      if(pdbfile.length()>0){
        for(unsigned i=0;i<pdb.size();++i){
          AtomNumber an=pdb.getAtomNumbers()[i];
          unsigned index=an.index();
          if( index>=unsigned(natoms) ) error("atom index in pdb exceeds the number of atoms in trajectory");
          masses[index]=pdb.getOccupancy()[i];
          charges[index]=pdb.getBeta()[i];
        }
      }
    } else if( checknatoms<0 && noatoms ){ 
      natoms=0; 
    }
    if( checknatoms<0 ){
      checknatoms=natoms;
      p.cmd("setNatoms",&natoms);
      p.cmd("init");
    }
    if(checknatoms!=natoms){
       std::string stepstr; Tools::convert(step,stepstr);
       error("number of atoms in frame " + stepstr + " does not match number of atoms in first frame");
    }

    coordinates.assign(3*natoms,real(0.0));
    forces.assign(3*natoms,real(0.0));
    cell.assign(9,real(0.0));
    virial.assign(9,real(0.0));

    if( first_step || rnd.U01()>0.5){
      if(debug_pd){
        int npe=intracomm.Get_size();
        vector<int> loc(npe,0);
        vector<int> start(npe,0);
        for(int i=0;i<npe-1;i++){
          int cc=(natoms*2*rnd.U01())/npe;
          if(start[i]+cc>natoms) cc=natoms-start[i];
          loc[i]=cc;
          start[i+1]=start[i]+loc[i];
        }
        loc[npe-1]=natoms-start[npe-1];
        intracomm.Bcast(loc,0);
        intracomm.Bcast(start,0);
        pd_nlocal=loc[intracomm.Get_rank()];
        pd_start=start[intracomm.Get_rank()];
        if(intracomm.Get_rank()==0){
          fprintf(out,"\nDRIVER: Reassigning particle decomposition\n");
          fprintf(out,"DRIVER: "); for(int i=0;i<npe;i++) fprintf(out,"%d ",loc[i]); printf("\n");
          fprintf(out,"DRIVER: "); for(int i=0;i<npe;i++) fprintf(out,"%d ",start[i]); printf("\n");
        }
        p.cmd("setAtomsNlocal",&pd_nlocal);
        p.cmd("setAtomsContiguous",&pd_start);
      } else if(debug_dd){
        int npe=intracomm.Get_size();
        int rank=intracomm.Get_rank();
        dd_charges.assign(natoms,0.0);
        dd_masses.assign(natoms,0.0);
        dd_gatindex.assign(natoms,-1);
        dd_g2l.assign(natoms,-1);
        dd_coordinates.assign(3*natoms,0.0);
        dd_forces.assign(3*natoms,0.0);
        dd_nlocal=0;
        for(int i=0;i<natoms;++i){
          double r=rnd.U01()*npe;
          int n; for(n=0;n<npe;n++) if(n+1>r)break;
          plumed_assert(n<npe);
          if(n==rank){
            dd_gatindex[dd_nlocal]=i;
            dd_g2l[i]=dd_nlocal;
            dd_charges[dd_nlocal]=charges[i];
            dd_masses[dd_nlocal]=masses[i];
            dd_nlocal++;
          }
        }
        if(intracomm.Get_rank()==0){
          fprintf(out,"\nDRIVER: Reassigning particle decomposition\n");
        }
        p.cmd("setAtomsNlocal",&dd_nlocal);
        p.cmd("setAtomsGatindex",&dd_gatindex[0]);
      }
    }

    int plumedStopCondition=0;
    p.cmd("setStep",&step);
    p.cmd("setStopFlag",&plumedStopCondition);
    if(!noatoms){
       if(trajectory_fmt=="xyz"){
         if(!Tools::getline(fp,line)) error("premature end of trajectory file");

         std::vector<double> celld(9,0.0);
         if(pbc_cli_given==false) {
           std::vector<std::string> words;
           words=Tools::getWords(line);
           if(words.size()==3){
             sscanf(line.c_str(),"%100lf %100lf %100lf",&celld[0],&celld[4],&celld[8]);
           } else if(words.size()==9){
             sscanf(line.c_str(),"%100lf %100lf %100lf %100lf %100lf %100lf %100lf %100lf %100lf",
                    &celld[0], &celld[1], &celld[2],
                    &celld[3], &celld[4], &celld[5],
                    &celld[6], &celld[7], &celld[8]);
           } else error("needed box in second line of xyz file");
         } else {			// from command line
           celld=pbc_cli_box;
         }
         for(unsigned i=0;i<9;i++)cell[i]=real(celld[i]);
       }
  	   int ddist=0;
       // Read coordinates
       for(int i=0;i<natoms;i++){
         bool ok=Tools::getline(fp,line);
         if(!ok) error("premature end of trajectory file");
         double cc[3];
         if(trajectory_fmt=="xyz"){
           char dummy[1000];
           std::sscanf(line.c_str(),"%999s %100lf %100lf %100lf",dummy,&cc[0],&cc[1],&cc[2]);
         } else if(trajectory_fmt=="gro"){
           // do the gromacs way
           if(!i){
        	   //
        	   // calculate the distance between dots (as in gromacs gmxlib/confio.c, routine get_w_conf )
        	   //
        	   const char      *p1, *p2, *p3;
        	   p1 = strchr(line.c_str(), '.');
        	   if (p1 == NULL) error("seems there are no coordinates in the gro file");
        	   p2 = strchr(&p1[1], '.');
        	   if (p2 == NULL) error("seems there is only one coordinates in the gro file");
        	   ddist = p2 - p1;
        	   p3 = strchr(&p2[1], '.');
        	   if (p3 == NULL)error("seems there are only two coordinates in the gro file");
        	   if (p3 - p2 != ddist)error("not uniform spacing in fields in the gro file");
           }
           Tools::convert(line.substr(20,ddist),cc[0]);
           Tools::convert(line.substr(20+ddist,ddist),cc[1]);
           Tools::convert(line.substr(20+ddist+ddist,ddist),cc[2]);
         } else plumed_error();
         if(!debug_pd || ( i>=pd_start && i<pd_start+pd_nlocal) ){
           coordinates[3*i]=real(cc[0]);
           coordinates[3*i+1]=real(cc[1]);
           coordinates[3*i+2]=real(cc[2]);
         }
       }
       if(trajectory_fmt=="gro"){
         if(!Tools::getline(fp,line)) error("premature end of trajectory file");
         std::vector<string> words=Tools::getWords(line);
         if(words.size()<3) error("cannot understand box format");
         Tools::convert(words[0],cell[0]);
         Tools::convert(words[1],cell[4]);
         Tools::convert(words[2],cell[8]);
         if(words.size()>3) Tools::convert(words[3],cell[1]);
         if(words.size()>4) Tools::convert(words[4],cell[2]);
         if(words.size()>5) Tools::convert(words[5],cell[3]);
         if(words.size()>6) Tools::convert(words[6],cell[5]);
         if(words.size()>7) Tools::convert(words[7],cell[6]);
         if(words.size()>8) Tools::convert(words[8],cell[7]);
       }

       if(debug_dd){
         for(int i=0;i<dd_nlocal;++i){
           int kk=dd_gatindex[i];
           dd_coordinates[3*i+0]=coordinates[3*kk+0];
           dd_coordinates[3*i+1]=coordinates[3*kk+1];
           dd_coordinates[3*i+2]=coordinates[3*kk+2];
         }
         p.cmd("setForces",&dd_forces[0]);
         p.cmd("setPositions",&dd_coordinates[0]);
         p.cmd("setMasses",&dd_masses[0]);
         p.cmd("setCharges",&dd_charges[0]);
       } else {
         p.cmd("setForces",&forces[3*pd_start]);
         p.cmd("setPositions",&coordinates[3*pd_start]);
         p.cmd("setMasses",&masses[pd_start]);
         p.cmd("setCharges",&charges[pd_start]);
       }
       p.cmd("setBox",&cell[0]);
       p.cmd("setVirial",&virial[0]);
   }
   p.cmd("calc");

// this is necessary as only processor zero is adding to the virial:
   intracomm.Bcast(virial,0);
   if(debug_pd) intracomm.Sum(forces);
   if(debug_dd){
     for(int i=0;i<dd_nlocal;i++){
       forces[3*dd_gatindex[i]+0]=dd_forces[3*i+0];
       forces[3*dd_gatindex[i]+1]=dd_forces[3*i+1];
       forces[3*dd_gatindex[i]+2]=dd_forces[3*i+2];
     }
     dd_forces.assign(3*natoms,0.0);
     intracomm.Sum(forces);
   }
   if(debug_grex &&step%grex_stride==0){
     p.cmd("GREX savePositions");
     if(intracomm.Get_rank()>0){
       p.cmd("GREX prepare");
     } else {
       int r=intercomm.Get_rank();
       int n=intercomm.Get_size();
       int partner=r+(2*((r+step/grex_stride)%2))-1;
       if(partner<0)partner=0;
       if(partner>=n) partner=n-1;
       p.cmd("GREX setPartner",&partner);
       p.cmd("GREX calculate");
       p.cmd("GREX shareAllDeltaBias");
       for(int i=0;i<n;i++){
         string s; Tools::convert(i,s);
         real a; s="GREX getDeltaBias "+s; p.cmd(s.c_str(),&a);
         if(grex_log) fprintf(grex_log," %f",a);
       }
       if(grex_log) fprintf(grex_log,"\n");
     }
   }


   if(fp_forces){
     fprintf(fp_forces,"%d\n",natoms);
     string fmt=dumpforcesFmt+" "+dumpforcesFmt+" "+dumpforcesFmt+"\n";
     fprintf(fp_forces,fmt.c_str(),virial[0],virial[4],virial[8]);
     fmt="X "+fmt;
     for(int i=0;i<natoms;i++)
       fprintf(fp_forces,fmt.c_str(),forces[3*i],forces[3*i+1],forces[3*i+2]);
   }

    if(noatoms && plumedStopCondition) break;

    step+=stride;
  }
  p.cmd("runFinalJobs");

  if(fp_forces) fclose(fp_forces);
  if(fp && fp!=in)fclose(fp);
  if(grex_log) fclose(grex_log);
  
  return 0;
}