Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2011-2023 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 : namespace PLMD {
31 :
32 208 : GREX::GREX(PlumedMain&p):
33 208 : initialized(false),
34 208 : plumedMain(p),
35 208 : atoms(p.getAtoms()),
36 208 : partner(-1), // = unset
37 208 : localDeltaBias(0),
38 208 : foreignDeltaBias(0),
39 208 : localUNow(0),
40 208 : localUSwap(0),
41 208 : myreplica(-1) { // = unset
42 208 : p.setSuffix(".NA");
43 208 : }
44 :
45 416 : GREX::~GREX() {
46 : // empty destructor to delete unique_ptr
47 416 : }
48 :
49 : #define CHECK_INIT(ini,word) plumed_massert(ini,"cmd(\"" + word +"\") should be only used after GREX initialization")
50 : #define CHECK_NOTINIT(ini,word) plumed_massert(!(ini),"cmd(\"" + word +"\") should be only used before GREX initialization")
51 : #define CHECK_NOTNULL(val,word) plumed_massert(val,"NULL pointer received in cmd(\"GREX " + word + "\")");
52 :
53 1089 : void GREX::cmd(const std::string&key,const TypesafePtr & val) {
54 : // Enumerate all possible commands:
55 : enum {
56 : #include "GREXEnum.inc"
57 : };
58 :
59 : // Static object (initialized once) containing the map of commands:
60 : const static std::unordered_map<std::string, int> word_map = {
61 : #include "GREXMap.inc"
62 4574 : };
63 :
64 1089 : std::vector<std::string> words=Tools::getWords(key);
65 1089 : unsigned nw=words.size();
66 1089 : if(nw==0) {
67 : // do nothing
68 : } else {
69 : int iword=-1;
70 : const auto it=word_map.find(words[0]);
71 1089 : if(it!=word_map.end()) {
72 1089 : iword=it->second;
73 : }
74 1089 : switch(iword) {
75 : case cmd_initialized:
76 0 : CHECK_NOTNULL(val,key);
77 0 : val.set(int(initialized));
78 : break;
79 208 : case cmd_setMPIIntracomm:
80 208 : CHECK_NOTINIT(initialized,key);
81 416 : intracomm.Set_comm(val.get<const void*>());
82 208 : break;
83 160 : case cmd_setMPIIntercomm:
84 160 : CHECK_NOTINIT(initialized,key);
85 320 : intercomm.Set_comm(val.get<const void*>());
86 320 : plumedMain.multi_sim_comm.Set_comm(val.get<const void*>());
87 160 : break;
88 0 : case cmd_setMPIFIntracomm:
89 0 : CHECK_NOTINIT(initialized,key);
90 0 : intracomm.Set_fcomm(val.get<const void*>());
91 0 : break;
92 0 : case cmd_setMPIFIntercomm:
93 0 : CHECK_NOTINIT(initialized,key);
94 0 : intercomm.Set_fcomm(val.get<const void*>());
95 0 : plumedMain.multi_sim_comm.Set_fcomm(val.get<const void*>());
96 0 : break;
97 208 : case cmd_init:
98 208 : CHECK_NOTINIT(initialized,key);
99 208 : initialized=true;
100 : // note that for PEs!=root this is automatically 0 (comm defaults to MPI_COMM_SELF)
101 208 : myreplica=intercomm.Get_rank();
102 208 : intracomm.Sum(myreplica);
103 : {
104 : std::string s;
105 208 : Tools::convert(myreplica,s);
106 416 : plumedMain.setSuffix("."+s);
107 : }
108 208 : break;
109 57 : case cmd_prepare:
110 57 : CHECK_INIT(initialized,key);
111 57 : if(intracomm.Get_rank()==0) {
112 : return;
113 : }
114 57 : intracomm.Bcast(partner,0);
115 57 : calculate();
116 : break;
117 57 : case cmd_setPartner:
118 57 : CHECK_INIT(initialized,key);
119 57 : partner=val.get<int>();
120 57 : break;
121 114 : case cmd_savePositions:
122 114 : CHECK_INIT(initialized,key);
123 114 : savePositions();
124 : break;
125 57 : case cmd_calculate:
126 57 : CHECK_INIT(initialized,key);
127 57 : if(intracomm.Get_rank()!=0) {
128 : return;
129 : }
130 57 : intracomm.Bcast(partner,0);
131 57 : calculate();
132 : break;
133 0 : case cmd_getLocalDeltaBias:
134 0 : CHECK_INIT(initialized,key);
135 0 : CHECK_NOTNULL(val,key);
136 0 : atoms.double2MD(localDeltaBias/(atoms.getMDUnits().getEnergy()/atoms.getUnits().getEnergy()),val);
137 0 : break;
138 0 : case cmd_cacheLocalUNow:
139 0 : CHECK_INIT(initialized,key);
140 0 : CHECK_NOTNULL(val,key);
141 : {
142 : double x;
143 0 : atoms.MD2double(val,x);
144 0 : localUNow=x*(atoms.getMDUnits().getEnergy()/atoms.getUnits().getEnergy());
145 0 : intracomm.Sum(localUNow);
146 : }
147 0 : break;
148 0 : case cmd_cacheLocalUSwap:
149 0 : CHECK_INIT(initialized,key);
150 0 : CHECK_NOTNULL(val,key);
151 : {
152 : double x;
153 0 : atoms.MD2double(val,x);
154 0 : localUSwap=x*(atoms.getMDUnits().getEnergy()/atoms.getUnits().getEnergy());
155 0 : intracomm.Sum(localUSwap);
156 : }
157 0 : break;
158 0 : case cmd_getForeignDeltaBias:
159 0 : CHECK_INIT(initialized,key);
160 0 : CHECK_NOTNULL(val,key);
161 0 : atoms.double2MD(foreignDeltaBias/(atoms.getMDUnits().getEnergy()/atoms.getUnits().getEnergy()),val);
162 0 : break;
163 57 : case cmd_shareAllDeltaBias:
164 57 : CHECK_INIT(initialized,key);
165 57 : if(intracomm.Get_rank()!=0) {
166 : return;
167 : }
168 57 : allDeltaBias.assign(intercomm.Get_size(),0.0);
169 57 : allDeltaBias[intercomm.Get_rank()]=localDeltaBias;
170 57 : intercomm.Sum(allDeltaBias);
171 : break;
172 171 : case cmd_getDeltaBias:
173 171 : CHECK_INIT(initialized,key);
174 171 : CHECK_NOTNULL(val,key);
175 171 : plumed_assert(nw==2);
176 171 : plumed_massert(allDeltaBias.size()==static_cast<unsigned>(intercomm.Get_size()),
177 : "to retrieve bias with cmd(\"GREX getDeltaBias\"), first share it with cmd(\"GREX shareAllDeltaBias\")");
178 : {
179 : unsigned rep;
180 171 : Tools::convert(words[1],rep);
181 171 : plumed_massert(rep<allDeltaBias.size(),"replica index passed to cmd(\"GREX getDeltaBias\") is out of range");
182 171 : double d=allDeltaBias[rep]/(atoms.getMDUnits().getEnergy()/atoms.getUnits().getEnergy());
183 171 : atoms.double2MD(d,val);
184 : }
185 171 : break;
186 0 : default:
187 0 : plumed_merror("cannot interpret cmd(\" GREX" + key + "\"). check plumed developers manual to see the available commands.");
188 : break;
189 : }
190 : }
191 1089 : }
192 :
193 114 : void GREX::savePositions() {
194 114 : plumedMain.prepareDependencies();
195 114 : plumedMain.resetActive(true);
196 114 : atoms.shareAll();
197 114 : plumedMain.waitData();
198 114 : std::ostringstream o;
199 114 : atoms.writeBinary(o);
200 114 : buffer=o.str();
201 114 : }
202 :
203 114 : void GREX::calculate() {
204 : unsigned nn=buffer.size();
205 114 : std::vector<char> rbuf(nn);
206 114 : localDeltaBias=-plumedMain.getBias();
207 114 : if(intracomm.Get_rank()==0) {
208 57 : Communicator::Request req=intercomm.Isend(buffer,partner,1066);
209 57 : intercomm.Recv(rbuf,partner,1066);
210 57 : req.wait();
211 : }
212 114 : intracomm.Bcast(rbuf,0);
213 114 : std::istringstream i(std::string(&rbuf[0],rbuf.size()));
214 114 : atoms.readBinary(i);
215 114 : plumedMain.setExchangeStep(true);
216 114 : plumedMain.prepareDependencies();
217 114 : plumedMain.justCalculate();
218 114 : plumedMain.setExchangeStep(false);
219 114 : localDeltaBias+=plumedMain.getBias();
220 114 : localDeltaBias+=localUSwap-localUNow;
221 114 : if(intracomm.Get_rank()==0) {
222 57 : Communicator::Request req=intercomm.Isend(localDeltaBias,partner,1067);
223 57 : intercomm.Recv(foreignDeltaBias,partner,1067);
224 57 : req.wait();
225 : }
226 114 : intracomm.Bcast(foreignDeltaBias,0);
227 228 : }
228 :
229 : }
|