Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2011-2019 The plumed team
3 : (see the PEOPLE file at the root of the distribution for a list of names)
4 :
5 : See http://www.plumed.org for more information.
6 :
7 : This file is part of plumed, version 2.
8 :
9 : plumed is free software: you can redistribute it and/or modify
10 : it under the terms of the GNU Lesser General Public License as published by
11 : the Free Software Foundation, either version 3 of the License, or
12 : (at your option) any later version.
13 :
14 : plumed is distributed in the hope that it will be useful,
15 : but WITHOUT ANY WARRANTY; without even the implied warranty of
16 : MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 : GNU Lesser General Public License for more details.
18 :
19 : You should have received a copy of the GNU Lesser General Public License
20 : along with plumed. If not, see <http://www.gnu.org/licenses/>.
21 : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
22 : #include "GREX.h"
23 : #include "PlumedMain.h"
24 : #include "Atoms.h"
25 : #include "tools/Tools.h"
26 : #include "tools/Communicator.h"
27 : #include <sstream>
28 : #include <unordered_map>
29 :
30 : using namespace std;
31 : namespace PLMD {
32 :
33 186 : GREX::GREX(PlumedMain&p):
34 : initialized(false),
35 : plumedMain(p),
36 : atoms(p.getAtoms()),
37 : partner(-1), // = unset
38 : localDeltaBias(0),
39 : foreignDeltaBias(0),
40 : localUNow(0),
41 : localUSwap(0),
42 930 : myreplica(-1) // = unset
43 : {
44 372 : p.setSuffix(".NA");
45 186 : }
46 :
47 372 : GREX::~GREX() {
48 : // empty destructor to delete unique_ptr
49 372 : }
50 :
51 : #define CHECK_INIT(ini,word) plumed_massert(ini,"cmd(\"" + word +"\") should be only used after GREX initialization")
52 : #define CHECK_NOTINIT(ini,word) plumed_massert(!(ini),"cmd(\"" + word +"\") should be only used before GREX initialization")
53 : #define CHECK_NOTNULL(val,word) plumed_massert(val,"NULL pointer received in cmd(\"GREX " + word + "\")");
54 :
55 1023 : void GREX::cmd(const string&key,void*val) {
56 :
57 : // Enumerate all possible commands:
58 : enum {
59 : #include "GREXEnum.inc"
60 : };
61 :
62 : // Static object (initialized once) containing the map of commands:
63 : const static std::unordered_map<std::string, int> word_map = {
64 : #include "GREXMap.inc"
65 1389 : };
66 :
67 1023 : std::vector<std::string> words=Tools::getWords(key);
68 1023 : unsigned nw=words.size();
69 1023 : if(nw==0) {
70 : // do nothing
71 : } else {
72 : int iword=-1;
73 : const auto it=word_map.find(words[0]);
74 1023 : if(it!=word_map.end()) iword=it->second;
75 1023 : switch(iword) {
76 : case cmd_initialized:
77 0 : CHECK_NOTNULL(val,key);
78 0 : *static_cast<int*>(val)=initialized;
79 0 : break;
80 : case cmd_setMPIIntracomm:
81 186 : CHECK_NOTINIT(initialized,key);
82 186 : intracomm.Set_comm(val);
83 : break;
84 : case cmd_setMPIIntercomm:
85 138 : CHECK_NOTINIT(initialized,key);
86 138 : intercomm.Set_comm(val);
87 138 : plumedMain.multi_sim_comm.Set_comm(val);
88 : break;
89 : case cmd_setMPIFIntracomm:
90 0 : CHECK_NOTINIT(initialized,key);
91 0 : intracomm.Set_fcomm(val);
92 : break;
93 : case cmd_setMPIFIntercomm:
94 0 : CHECK_NOTINIT(initialized,key);
95 0 : intercomm.Set_fcomm(val);
96 0 : plumedMain.multi_sim_comm.Set_fcomm(val);
97 : break;
98 : case cmd_init:
99 186 : CHECK_NOTINIT(initialized,key);
100 186 : initialized=true;
101 : // note that for PEs!=root this is automatically 0 (comm defaults to MPI_COMM_SELF)
102 186 : myreplica=intercomm.Get_rank();
103 186 : intracomm.Sum(myreplica);
104 : {
105 : std::string s;
106 186 : Tools::convert(myreplica,s);
107 372 : plumedMain.setSuffix("."+s);
108 : }
109 186 : break;
110 : case cmd_prepare:
111 57 : CHECK_INIT(initialized,key);
112 57 : if(intracomm.Get_rank()==0) return;
113 57 : intracomm.Bcast(partner,0);
114 57 : calculate();
115 : break;
116 : case cmd_setPartner:
117 57 : CHECK_INIT(initialized,key);
118 57 : partner=*static_cast<int*>(val);
119 57 : break;
120 : case cmd_savePositions:
121 114 : CHECK_INIT(initialized,key);
122 114 : savePositions();
123 : break;
124 : case cmd_calculate:
125 57 : CHECK_INIT(initialized,key);
126 57 : if(intracomm.Get_rank()!=0) return;
127 57 : intracomm.Bcast(partner,0);
128 57 : calculate();
129 : break;
130 : case cmd_getLocalDeltaBias:
131 0 : CHECK_INIT(initialized,key);
132 0 : CHECK_NOTNULL(val,key);
133 0 : atoms.double2MD(localDeltaBias/(atoms.getMDUnits().getEnergy()/atoms.getUnits().getEnergy()),val);
134 0 : break;
135 : case cmd_cacheLocalUNow:
136 0 : CHECK_INIT(initialized,key);
137 0 : CHECK_NOTNULL(val,key);
138 : {
139 : double x;
140 0 : atoms.MD2double(val,x);
141 0 : localUNow=x*(atoms.getMDUnits().getEnergy()/atoms.getUnits().getEnergy());
142 0 : intracomm.Sum(localUNow);
143 : }
144 0 : break;
145 : case cmd_cacheLocalUSwap:
146 0 : CHECK_INIT(initialized,key);
147 0 : CHECK_NOTNULL(val,key);
148 : {
149 : double x;
150 0 : atoms.MD2double(val,x);
151 0 : localUSwap=x*(atoms.getMDUnits().getEnergy()/atoms.getUnits().getEnergy());
152 0 : intracomm.Sum(localUSwap);
153 : }
154 0 : break;
155 : case cmd_getForeignDeltaBias:
156 0 : CHECK_INIT(initialized,key);
157 0 : CHECK_NOTNULL(val,key);
158 0 : atoms.double2MD(foreignDeltaBias/(atoms.getMDUnits().getEnergy()/atoms.getUnits().getEnergy()),val);
159 0 : break;
160 : case cmd_shareAllDeltaBias:
161 57 : CHECK_INIT(initialized,key);
162 57 : if(intracomm.Get_rank()!=0) return;
163 114 : allDeltaBias.assign(intercomm.Get_size(),0.0);
164 114 : allDeltaBias[intercomm.Get_rank()]=localDeltaBias;
165 57 : intercomm.Sum(allDeltaBias);
166 : break;
167 : case cmd_getDeltaBias:
168 171 : CHECK_INIT(initialized,key);
169 171 : CHECK_NOTNULL(val,key);
170 171 : plumed_assert(nw==2);
171 171 : plumed_massert(allDeltaBias.size()==static_cast<unsigned>(intercomm.Get_size()),
172 0 : "to retrieve bias with cmd(\"GREX getDeltaBias\"), first share it with cmd(\"GREX shareAllDeltaBias\")");
173 : {
174 : unsigned rep;
175 171 : Tools::convert(words[1],rep);
176 342 : plumed_massert(rep<allDeltaBias.size(),"replica index passed to cmd(\"GREX getDeltaBias\") is out of range");
177 171 : double d=allDeltaBias[rep]/(atoms.getMDUnits().getEnergy()/atoms.getUnits().getEnergy());
178 171 : atoms.double2MD(d,val);
179 : }
180 171 : break;
181 : default:
182 0 : plumed_merror("cannot interpret cmd(\" GREX" + key + "\"). check plumed developers manual to see the available commands.");
183 : break;
184 : }
185 1023 : }
186 : }
187 :
188 114 : void GREX::savePositions() {
189 114 : plumedMain.prepareDependencies();
190 114 : plumedMain.resetActive(true);
191 114 : atoms.shareAll();
192 114 : plumedMain.waitData();
193 114 : ostringstream o;
194 114 : atoms.writeBinary(o);
195 228 : buffer=o.str();
196 114 : }
197 :
198 114 : void GREX::calculate() {
199 : //fprintf(stderr,"CALCULATE %d %d\n",intercomm.Get_rank(),partner);
200 : unsigned nn=buffer.size();
201 114 : vector<char> rbuf(nn);
202 114 : localDeltaBias=-plumedMain.getBias();
203 114 : if(intracomm.Get_rank()==0) {
204 57 : Communicator::Request req=intercomm.Isend(buffer,partner,1066);
205 57 : intercomm.Recv(rbuf,partner,1066);
206 57 : req.wait();
207 : }
208 114 : intracomm.Bcast(rbuf,0);
209 342 : istringstream i(string(&rbuf[0],rbuf.size()));
210 114 : atoms.readBinary(i);
211 114 : plumedMain.setExchangeStep(true);
212 114 : plumedMain.prepareDependencies();
213 114 : plumedMain.justCalculate();
214 114 : plumedMain.setExchangeStep(false);
215 114 : localDeltaBias+=plumedMain.getBias();
216 114 : localDeltaBias+=localUSwap-localUNow;
217 114 : if(intracomm.Get_rank()==0) {
218 57 : Communicator::Request req=intercomm.Isend(localDeltaBias,partner,1067);
219 57 : intercomm.Recv(foreignDeltaBias,partner,1067);
220 57 : req.wait();
221 : //fprintf(stderr,">>> %d %d %20.12f %20.12f %20.12f %20.12f\n",intercomm.Get_rank(),partner,localDeltaBias,foreignDeltaBias,localUSwap,localUNow);
222 : }
223 114 : intracomm.Bcast(foreignDeltaBias,0);
224 114 : }
225 :
226 5874 : }
|