Skip to content

Commit 127f9e2

Browse files
committed
New tool poisson_integrator_conv
1 parent 1ee3ec1 commit 127f9e2

File tree

6 files changed

+154
-14
lines changed

6 files changed

+154
-14
lines changed

colvartools/Makefile

+13-3
Original file line numberDiff line numberDiff line change
@@ -2,26 +2,36 @@ SHELL=/bin/sh
22
#########################
33
# adjust as needed
44
# settings for linux with GCC
5+
COLVARS_SRC=../src
56
CXX=g++
6-
CXXFLAGS=-Wall -O2
7+
CXXFLAGS=-Wall -O2 -I$(COLVARS_SRC)
78
EXT=
89
# settings for mingw cross-compiler to windows (32-bit)
910
#CXX=i688-w64-mingw64-g++
1011
#CXXFLAGS=-Wall -O2
1112
#EXT=.exe
1213
#########################
1314

14-
all: abf_integrate$(EXT)
15+
all: poisson_integrator$(EXT) abf_integrate$(EXT)
1516

1617
clean:
17-
-rm *~ *.o abf_integrate$(EXT) *.exe
18+
-rm *~ *.o abf_integrate$(EXT) poisson_integrator$(EXT) *.exe
1819

1920
abf_integrate$(EXT): abf_integrate.o abf_data.o
2021
$(CXX) -o $@ $(CXXFLAGS) $^
2122

23+
poisson_integrator$(EXT): poisson_integrator.o libcolvars
24+
$(CXX) -o $@ $(CXXFLAGS) poisson_integrator.o $(COLVARS_SRC)/libcolvars.a
25+
26+
poisson_integrator_conv: poisson_integrator_conv.o libcolvars
27+
$(CXX) -o $@ $(CXXFLAGS) poisson_integrator_conv.o $(COLVARS_SRC)/libcolvars.a
28+
2229
%.o: %.cpp
2330
$(CXX) -o $@ -c $(CXXFLAGS) $<
2431

32+
libcolvars:
33+
$(MAKE) -C $(COLVARS_SRC) libcolvars.a
34+
2535
# dependencies
2636
abf_integrate.o: abf_integrate.cpp abf_data.h
2737
abf_data.o: abf_data.cpp abf_data.h

colvartools/README.md

+4-2
Original file line numberDiff line numberDiff line change
@@ -5,7 +5,9 @@ This directory contains both standalone tools and Colvars scripts.
55
## Standalone tools
66
| File name | Summary |
77
| ------------- | ------------- |
8-
| **abf_integrate** | Post-process gradient files produced by ABF and related methods, to generate a PMF. Superseded by builtin integration for dimensions 2 and 3, still needed for higher-dimension PMFs. Build using the provided **Makefile**.|
8+
| **abf_integrate** | Process free energy gradient from ABF/TI to generate a free energy surface using a legacy MCMC procedure. Superseded by Poisson integration (builtin or standalone) for dimensions 2 and 3, still needed for higher-dimension PMFs. Build using the provided **Makefile**.|
9+
| **poisson_integrator** | Process free energy gradient from ABF/TI to generate a free energy surface as the solution of a Poisson equation. Build using the provided **Makefile**.|
10+
| **poisson_integrator_conv** | Process free energy gradient from ABF/TI to generate a free energy surface as the solution of a Poisson equation, monitoring convergence towards a given discretized scalar field. Build using the provided **Makefile**.|
911
| **noe_to_colvars.py** | Parse an X-PLOR style list of assign commands for NOE restraints.|
1012
| **plot_colvars_traj.py** | Select variables from a Colvars trajectory file and optionally plot them as a 1D graph as a function of time or of one of the variables.|
1113
| **quaternion2rmatrix.tcl** | As the name says.|
@@ -16,4 +18,4 @@ This directory contains both standalone tools and Colvars scripts.
1618
| File name | Summary |
1719
| ------------- | ------------- |
1820
| **abmd.tcl** | Adiabatic Biased MD after Marchi & Ballone JCP 1999 / Paci & Karplus JMB 1999; implemented in 20 lines of Tcl.|
19-
| **pathCV.tcl** | Path-based collective variables after Branduardi et al. (2007). Optimized implementation that calculates only the necessary distances.|
21+
| **pathCV.tcl** | Path-based collective variables after Branduardi et al. (2007). Optimized implementation that calculates only the necessary distances.|

colvartools/poisson_integrator.cpp

+11-5
Original file line numberDiff line numberDiff line change
@@ -8,7 +8,7 @@
88
int main (int argc, char *argv[]) {
99

1010
if (argc < 2) {
11-
std::cerr << "One argument needed: file name.\n";
11+
std::cerr << "\n\nOne argument needed: gradient multicol file name.\n";
1212
return 1;
1313
}
1414

@@ -17,18 +17,24 @@ int main (int argc, char *argv[]) {
1717

1818
std::string gradfile (argv[1]);
1919
std::shared_ptr<colvar_grid_gradient> grad_ptr = std::make_shared<colvar_grid_gradient>(gradfile);
20+
if (cvm::get_error()) { return -1; }
2021

21-
int itmax = 1000;
22+
int itmax = 10000;
2223
cvm::real err;
23-
cvm::real tol = 1e-6;
24+
cvm::real tol = 1e-8;
2425

2526
integrate_potential potential(grad_ptr);
2627
potential.set_div();
2728
potential.integrate(itmax, tol, err);
2829
potential.set_zero_minimum();
2930

30-
potential.write_multicol(std::string(gradfile + ".int"),
31-
"integrated potential");
31+
if (potential.num_variables() < 3) {
32+
std::cout << "\nWriting integrated potential in multicol format to " + gradfile + ".int\n";
33+
potential.write_multicol(std::string(gradfile + ".int"), "integrated potential");
34+
} else { // Write 3D grids to more convenient DX format
35+
std::cout << "\nWriting integrated potential in OpenDX format to " + gradfile + ".int.dx\n";
36+
potential.write_opendx(std::string(gradfile + ".int.dx"), "integrated potential");
37+
}
3238

3339
delete colvars;
3440
return 0;
+64
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,64 @@
1+
#include <iostream>
2+
3+
#include "colvargrid.h"
4+
#include "colvarproxy.h"
5+
6+
// Integrate provided gradients while monitoring convergence towards a provided scalar grid
7+
// (typically the result of a previous integration)
8+
9+
int main (int argc, char *argv[])
10+
{
11+
colvarproxy *proxy = new colvarproxy();
12+
colvarmodule *colvars = new colvarmodule(proxy);
13+
14+
if (argc < 2) {
15+
std::cerr << "\n\nOne argument needed: gradient multicol file name.\n";
16+
return 1;
17+
}
18+
19+
std::string gradfile (argv[1]);
20+
std::shared_ptr<colvar_grid_gradient> grad_ptr = std::make_shared<colvar_grid_gradient>(gradfile);
21+
if (cvm::get_error()) { return -1; }
22+
23+
cvm::real err = 1.;
24+
cvm::real tol = 1e-10;
25+
26+
integrate_potential potential(grad_ptr);
27+
potential.set_div();
28+
29+
// Load reference
30+
colvar_grid_scalar ref(gradfile + ".ref");
31+
if (cvm::get_error()) { return -1; }
32+
33+
if (ref.number_of_points() != potential.number_of_points()) {
34+
cvm::error("Reference grid has wrong number of points: " + cvm::to_str(ref.number_of_points()) + "\n");
35+
return -1;
36+
}
37+
38+
// Monitor convergence
39+
int rounds = 100;
40+
int steps_per_round = 10;
41+
42+
std::cout << 0 << " " << 0 << " " << ref.grid_rmsd(potential) << std::endl;
43+
44+
for (int i = 0; i < rounds && err > tol; i++) {
45+
potential.reset(0.);
46+
potential.integrate(steps_per_round * (i+1), tol, err);
47+
potential.set_zero_minimum();
48+
49+
std::cout << (i+1)*steps_per_round << " " << err << " " << ref.grid_rmsd(potential) << std::endl;
50+
51+
char buff[100];
52+
snprintf(buff, sizeof(buff), "%04i", steps_per_round * (i+1));
53+
}
54+
55+
if (potential.num_variables() < 3) {
56+
std::cout << "\nWriting integrated potential in multicol format to " + gradfile + ".int\n";
57+
potential.write_multicol(std::string(gradfile + ".int"), "integrated potential");
58+
} else { // Write 3D grids to more convenient DX format
59+
std::cout << "\nWriting integrated potential in OpenDX format to " + gradfile + ".int.dx\n";
60+
potential.write_opendx(std::string(gradfile + ".int.dx"), "integrated potential");
61+
}
62+
delete colvars;
63+
return 0;
64+
}

src/colvargrid.cpp

+58-3
Original file line numberDiff line numberDiff line change
@@ -142,6 +142,62 @@ colvar_grid_scalar::colvar_grid_scalar(std::vector<colvar *> &colvars, bool marg
142142
{
143143
}
144144

145+
146+
colvar_grid_scalar::colvar_grid_scalar(std::string const &filename)
147+
: colvar_grid<cvm::real>(),
148+
samples(NULL)
149+
{
150+
std::istream &is = cvm::main()->proxy->input_stream(filename, "multicol scalar file");
151+
if (!is) {
152+
return;
153+
}
154+
155+
// Data in the header: nColvars, then for each
156+
// xiMin, dXi, nPoints, periodic flag
157+
158+
std::string hash;
159+
size_t i;
160+
161+
if ( !(is >> hash) || (hash != "#") ) {
162+
cvm::error("Error reading grid at position "+
163+
cvm::to_str(static_cast<size_t>(is.tellg()))+
164+
" in stream(read \"" + hash + "\")\n");
165+
return;
166+
}
167+
168+
is >> nd;
169+
mult = 1;
170+
std::vector<cvm::real> lower_in(nd), widths_in(nd);
171+
std::vector<int> nx_in(nd);
172+
std::vector<int> periodic_in(nd);
173+
174+
for (i = 0; i < nd; i++ ) {
175+
if ( !(is >> hash) || (hash != "#") ) {
176+
cvm::error("Error reading grid at position "+
177+
cvm::to_str(static_cast<size_t>(is.tellg()))+
178+
" in stream(read \"" + hash + "\")\n");
179+
return;
180+
}
181+
182+
is >> lower_in[i] >> widths_in[i] >> nx_in[i] >> periodic_in[i];
183+
}
184+
185+
this->setup(nx_in, 0., mult);
186+
187+
widths = widths_in;
188+
189+
for (i = 0; i < nd; i++ ) {
190+
lower_boundaries.push_back(colvarvalue(lower_in[i]));
191+
periodic.push_back(static_cast<bool>(periodic_in[i]));
192+
}
193+
194+
// Reset the istream for read_multicol, which expects the whole file
195+
is.clear();
196+
is.seekg(0);
197+
read_multicol(is);
198+
cvm::main()->proxy->close_input_stream(filename);
199+
}
200+
145201
colvar_grid_scalar::~colvar_grid_scalar()
146202
{
147203
}
@@ -351,12 +407,11 @@ colvar_grid_gradient::colvar_grid_gradient(std::vector<colvar *> &colvars, std::
351407
}
352408

353409

354-
colvar_grid_gradient::colvar_grid_gradient(std::string &filename)
410+
colvar_grid_gradient::colvar_grid_gradient(std::string const &filename)
355411
: colvar_grid<cvm::real>(),
356412
samples(NULL)
357413
{
358-
std::istream &is = cvm::main()->proxy->input_stream(filename,
359-
"gradient file");
414+
std::istream &is = cvm::main()->proxy->input_stream(filename, "gradient file");
360415
if (!is) {
361416
return;
362417
}

src/colvargrid.h

+4-1
Original file line numberDiff line numberDiff line change
@@ -1262,6 +1262,9 @@ class colvar_grid_scalar : public colvar_grid<cvm::real>
12621262
colvar_grid_scalar(std::vector<colvar *> &colvars,
12631263
bool add_extra_bin = false);
12641264

1265+
/// Constructor from a multicol file
1266+
colvar_grid_scalar(std::string const &filename);
1267+
12651268
/// Accumulate the value
12661269
inline void acc_value(std::vector<int> const &ix,
12671270
cvm::real const &new_value,
@@ -1573,7 +1576,7 @@ class colvar_grid_gradient : public colvar_grid<cvm::real>
15731576
colvar_grid_gradient(std::vector<colvar *> &colvars);
15741577

15751578
/// Constructor from a multicol file
1576-
colvar_grid_gradient(std::string &filename);
1579+
colvar_grid_gradient(std::string const &filename);
15771580

15781581
/// Constructor from a vector of colvars and a pointer to the count grid
15791582
colvar_grid_gradient(std::vector<colvar *> &colvars, std::shared_ptr<colvar_grid_count> samples_in);

0 commit comments

Comments
 (0)