Skip to content

Commit

Permalink
Removed resizing of vectors from MultiColvarTemplate::performTask
Browse files Browse the repository at this point in the history
This resizing is now done when MultiColvar is setup
  • Loading branch information
Gareth Aneurin Tribello authored and Gareth Aneurin Tribello committed Feb 3, 2025
1 parent c9df402 commit f3860b2
Show file tree
Hide file tree
Showing 7 changed files with 23 additions and 19 deletions.
2 changes: 1 addition & 1 deletion src/adjmat/TopologyMatrix.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -127,7 +127,7 @@ double TopologyMatrix::calculateWeight( const Vector& pos1, const Vector& pos2,
// Now run through all sea atoms
HistogramBead bead; bead.isNotPeriodic(); bead.setKernelType( kerneltype );
Vector g1derivf,g2derivf,lderivf; Tensor vir; double binlength = maxbins * binw_mat;
MultiValue tvals; tvals.resize( maxbins, myvals.getNumberOfDerivatives() );
MultiValue tvals; tvals.resize( maxbins, myvals.getNumberOfDerivatives(), 0 );
for(unsigned i=0; i<natoms; ++i) {
// Position of sea atom (this will be the origin)
Vector d2 = getPosition(i,myvals);
Expand Down
13 changes: 6 additions & 7 deletions src/colvar/MultiColvarTemplate.h
Original file line number Diff line number Diff line change
Expand Up @@ -42,6 +42,7 @@ class MultiColvarTemplate : public ActionWithVector {
static void registerKeywords(Keywords&);
explicit MultiColvarTemplate(const ActionOptions&);
unsigned getNumberOfDerivatives() override ;
unsigned getNumberOfAtomsPerTask() const override ;
void addValueWithDerivatives( const std::vector<unsigned>& shape=std::vector<unsigned>() ) override ;
void addComponentWithDerivatives( const std::string& name, const std::vector<unsigned>& shape=std::vector<unsigned>() ) override ;
void performTask( const unsigned&, MultiValue& ) const override ;
Expand Down Expand Up @@ -127,11 +128,15 @@ void MultiColvarTemplate<T>::addComponentWithDerivatives( const std::string& nam
std::vector<unsigned> s(1); s[0]=ablocks[0].size(); addComponent( name, s );
}

template <class T>
unsigned MultiColvarTemplate<T>::getNumberOfAtomsPerTask() const {
return ablocks.size();
}

template <class T>
void MultiColvarTemplate<T>::performTask( const unsigned& task_index, MultiValue& myvals ) const {
// Retrieve the positions
std::vector<Vector> & fpositions( myvals.getFirstAtomVector() );
if( fpositions.size()!=ablocks.size() ) fpositions.resize( ablocks.size() );
for(unsigned i=0; i<ablocks.size(); ++i) fpositions[i] = getPosition( ablocks[i][task_index] );
// If we are using pbc make whole
if( usepbc ) {
Expand All @@ -145,19 +150,13 @@ void MultiColvarTemplate<T>::performTask( const unsigned& task_index, MultiValue
}
} else if( fpositions.size()==1 ) fpositions[0]=delta(Vector(0.0,0.0,0.0),getPosition( ablocks[0][task_index] ) );
// Retrieve the masses and charges
myvals.resizeTemporyVector(2);
std::vector<double> & mass( myvals.getTemporyVector(0) );
std::vector<double> & charge( myvals.getTemporyVector(1) );
if( mass.size()!=ablocks.size() ) { mass.resize(ablocks.size()); charge.resize(ablocks.size()); }
for(unsigned i=0; i<ablocks.size(); ++i) { mass[i]=getMass( ablocks[i][task_index] ); charge[i]=getCharge( ablocks[i][task_index] ); }
// Make some space to store various things
std::vector<double> values( getNumberOfComponents() );
std::vector<Tensor> & virial( myvals.getFirstAtomVirialVector() );
std::vector<std::vector<Vector> > & derivs( myvals.getFirstAtomDerivativeVector() );
if( derivs.size()!=values.size() ) { derivs.resize( values.size() ); virial.resize( values.size() ); }
for(unsigned i=0; i<derivs.size(); ++i) {
if( derivs[i].size()<ablocks.size() ) derivs[i].resize( ablocks.size() );
}
// Calculate the CVs using the method in the Colvar
T::calculateCV( mode, mass, charge, fpositions, values, derivs, virial, this );
for(unsigned i=0; i<values.size(); ++i) myvals.setValue( i, values[i] );
Expand Down
8 changes: 6 additions & 2 deletions src/core/ActionWithVector.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -222,7 +222,9 @@ void ActionWithVector::runAllTasks() {
std::vector<double> omp_buffer;
const unsigned t=OpenMP::getThreadNum();
if( nt>1 ) omp_buffer.resize( bufsize, 0.0 );
if( myvals[t].getNumberOfValues()!=getNumberOfComponents() || myvals[t].getNumberOfDerivatives()!=nderivatives ) myvals[t].resize( getNumberOfComponents(), nderivatives );
if( myvals[t].getNumberOfValues()!=getNumberOfComponents() || myvals[t].getNumberOfDerivatives()!=nderivatives || myvals[t].getAtomVector().size()!=getNumberOfAtomsPerTask() ) {
myvals[t].resize( getNumberOfComponents(), nderivatives, getNumberOfAtomsPerTask() );
}
myvals[t].clearAll();

#pragma omp for nowait
Expand Down Expand Up @@ -379,7 +381,9 @@ bool ActionWithVector::checkForForces() {
if( omp_forces[t].size()!=forcesForApply.size() ) omp_forces[t].resize( forcesForApply.size(), 0.0 );
else omp_forces[t].assign( forcesForApply.size(), 0.0 );
}
if( myvals[t].getNumberOfValues()!=getNumberOfComponents() || myvals[t].getNumberOfDerivatives()!=nderiv ) myvals[t].resize( getNumberOfComponents(), nderiv );
if( myvals[t].getNumberOfValues()!=getNumberOfComponents() || myvals[t].getNumberOfDerivatives()!=nderiv || myvals[t].getAtomVector().size()!=getNumberOfAtomsPerTask() ) {
myvals[t].resize( getNumberOfComponents(), nderiv, getNumberOfAtomsPerTask() );
}
myvals[t].clearAll();

#pragma omp for nowait
Expand Down
2 changes: 2 additions & 0 deletions src/core/ActionWithVector.h
Original file line number Diff line number Diff line change
Expand Up @@ -84,6 +84,8 @@ class ActionWithVector:
/// Check if a mask has been set
int getNumberOfMasks() const ;
void calculateNumericalDerivatives(ActionWithValue* av) override;
/// This is for resizing the task list
virtual unsigned getNumberOfAtomsPerTask() const { return 0; }
/// Turn off the calculation of the derivatives during the forward pass through a calculation
bool doNotCalculateDerivatives() const override ;
/// Get the list of tasks that are active
Expand Down
2 changes: 1 addition & 1 deletion src/refdist/MatrixProductDiagonal.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -127,7 +127,7 @@ void MatrixProductDiagonal::performTask( const unsigned& task_index, MultiValue&
void MatrixProductDiagonal::calculate() {
if( getPntrToArgument(1)->getRank()==1 ) {
unsigned nder = getNumberOfDerivatives();
MultiValue myvals; myvals.resize( 1, nder ); performTask( 0, myvals );
MultiValue myvals; myvals.resize( 1, nder, 0 ); performTask( 0, myvals );

Value* myval=getPntrToComponent(0); myval->set( myvals.get(0) );
for(unsigned i=0; i<nder; ++i) myval->setDerivative( i, myvals.getDerivative(0,i) );
Expand Down
7 changes: 6 additions & 1 deletion src/tools/MultiValue.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -24,12 +24,17 @@

namespace PLMD {

void MultiValue::resize( const size_t& nvals, const size_t& nder ) {
void MultiValue::resize( const size_t& nvals, const size_t& nder, const size_t& natoms ) {
if( values.size()==nvals && nderivatives>nder ) return;
values.resize(nvals); nderivatives=nder; derivatives.resize( nvals*nder );
hasderiv.resize(nvals*nder,false); nactive.resize(nvals); active_list.resize(nvals*nder);
matrix_force_stash.resize(nder,0);
matrix_row_nderivatives=0; matrix_row_derivative_indices.resize(nder); atLeastOneSet=false;
if( natoms>0 ) {
tmp_vectors.resize(2); tmp_vectors[0].resize(natoms); tmp_vectors[1].resize(natoms);
tmp_atom_der.resize(nvals); tmp_atom_virial.resize(nvals); tmp_atoms.resize(natoms);
for(unsigned i=0; i<nvals; ++i) tmp_atom_der[i].resize(natoms);
}
}

void MultiValue::clearAll() {
Expand Down
8 changes: 1 addition & 7 deletions src/tools/MultiValue.h
Original file line number Diff line number Diff line change
Expand Up @@ -63,7 +63,7 @@ class MultiValue {
std::vector<std::vector<double> > tmp_vectors;
public:
MultiValue() : task_index(0), task2_index(0), nderivatives(0), atLeastOneSet(false), vector_call(false), nindices(0), nsplit(0), matrix_row_nderivatives(0) {}
void resize( const std::size_t& nvals, const std::size_t& nder );
void resize( const std::size_t& nvals, const std::size_t& nder, const std::size_t& natoms );
/// Set the task index prior to the loop
void setTaskIndex( const std::size_t& tindex );
///
Expand Down Expand Up @@ -91,7 +91,6 @@ class MultiValue {
std::vector<std::vector<Vector> >& getFirstAtomDerivativeVector();
const std::vector<std::vector<Vector> >& getConstFirstAtomDerivativeVector() const ;
std::vector<Tensor>& getFirstAtomVirialVector();
void resizeTemporyVector(const unsigned& n );
std::vector<double>& getTemporyVector(const unsigned& ind );
///
bool inVectorCall() const ;
Expand Down Expand Up @@ -322,11 +321,6 @@ double MultiValue::getStashedMatrixForce( const unsigned& jind ) const {
return matrix_force_stash[jind];
}

inline
void MultiValue::resizeTemporyVector(const unsigned& n ) {
if( n>tmp_vectors.size() ) tmp_vectors.resize(n);
}

inline
std::vector<double>& MultiValue::getTemporyVector(const unsigned& ind ) {
plumed_dbg_assert( ind<tmp_vectors.size() );
Expand Down

1 comment on commit f3860b2

@PlumedBot
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Found broken examples in automatic/ANGLES.tmp
Found broken examples in automatic/ANN.tmp
Found broken examples in automatic/CAVITY.tmp
Found broken examples in automatic/CLASSICAL_MDS.tmp
Found broken examples in automatic/CLUSTER_DIAMETER.tmp
Found broken examples in automatic/CLUSTER_DISTRIBUTION.tmp
Found broken examples in automatic/CLUSTER_PROPERTIES.tmp
Found broken examples in automatic/CONSTANT.tmp
Found broken examples in automatic/CONTACT_MATRIX.tmp
Found broken examples in automatic/CONTACT_MATRIX_PROPER.tmp
Found broken examples in automatic/COORDINATIONNUMBER.tmp
Found broken examples in automatic/DFSCLUSTERING.tmp
Found broken examples in automatic/DISTANCE_FROM_CONTOUR.tmp
Found broken examples in automatic/EDS.tmp
Found broken examples in automatic/EMMI.tmp
Found broken examples in automatic/ENVIRONMENTSIMILARITY.tmp
Found broken examples in automatic/FIND_CONTOUR.tmp
Found broken examples in automatic/FIND_CONTOUR_SURFACE.tmp
Found broken examples in automatic/FIND_SPHERICAL_CONTOUR.tmp
Found broken examples in automatic/FOURIER_TRANSFORM.tmp
Found broken examples in automatic/FUNCPATHGENERAL.tmp
Found broken examples in automatic/FUNCPATHMSD.tmp
Found broken examples in automatic/FUNNEL.tmp
Found broken examples in automatic/FUNNEL_PS.tmp
Found broken examples in automatic/GHBFIX.tmp
Found broken examples in automatic/GPROPERTYMAP.tmp
Found broken examples in automatic/HBOND_MATRIX.tmp
Found broken examples in automatic/INCLUDE.tmp
Found broken examples in automatic/INCYLINDER.tmp
Found broken examples in automatic/INENVELOPE.tmp
Found broken examples in automatic/INTERPOLATE_GRID.tmp
Found broken examples in automatic/LOCAL_AVERAGE.tmp
Found broken examples in automatic/MAZE_OPTIMIZER_BIAS.tmp
Found broken examples in automatic/MAZE_RANDOM_ACCELERATION_MD.tmp
Found broken examples in automatic/MAZE_SIMULATED_ANNEALING.tmp
Found broken examples in automatic/MAZE_STEERED_MD.tmp
Found broken examples in automatic/METATENSOR.tmp
Found broken examples in automatic/MULTICOLVARDENS.tmp
Found broken examples in automatic/OUTPUT_CLUSTER.tmp
Found broken examples in automatic/PAMM.tmp
Found broken examples in automatic/PCA.tmp
Found broken examples in automatic/PCAVARS.tmp
Found broken examples in automatic/PIV.tmp
Found broken examples in automatic/PLUMED.tmp
Found broken examples in automatic/PYCVINTERFACE.tmp
Found broken examples in automatic/PYTHONFUNCTION.tmp
Found broken examples in automatic/Q3.tmp
Found broken examples in automatic/Q4.tmp
Found broken examples in automatic/Q6.tmp
Found broken examples in automatic/QUATERNION.tmp
Found broken examples in automatic/SIZESHAPE_POSITION_LINEAR_PROJ.tmp
Found broken examples in automatic/SIZESHAPE_POSITION_MAHA_DIST.tmp
Found broken examples in automatic/SPRINT.tmp
Found broken examples in automatic/TETRAHEDRALPORE.tmp
Found broken examples in automatic/TORSIONS.tmp
Found broken examples in automatic/WHAM_WEIGHTS.tmp
Found broken examples in AnalysisPP.md
Found broken examples in CollectiveVariablesPP.md
Found broken examples in MiscelaneousPP.md

Please sign in to comment.