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