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 "tools/Tools.h"
25 : #include "tools/Communicator.h"
26 : #include <sstream>
27 :
28 : namespace PLMD {
29 :
30 208 : GREX::GREX(PlumedMain&p):
31 208 : initialized(false),
32 208 : plumedMain(p),
33 208 : partner(-1), // = unset
34 208 : localDeltaBias(0),
35 208 : foreignDeltaBias(0),
36 208 : localUNow(0),
37 208 : localUSwap(0),
38 208 : myreplica(-1) { // = unset
39 208 : p.setSuffix(".NA");
40 208 : }
41 :
42 416 : GREX::~GREX() {
43 : // empty destructor to delete unique_ptr
44 416 : }
45 :
46 : #define CHECK_INIT(ini,word) plumed_assert(ini) << "cmd(\"" << word <<"\") should be only used after GREX initialization"
47 : #define CHECK_NOTINIT(ini,word) plumed_assert(!(ini)) << "cmd(\"" << word <<"\") should be only used before GREX initialization"
48 : #define CHECK_NOTNULL(val,word) plumed_assert(val) << "NULL pointer received in cmd(\"GREX " << word << "\")"
49 :
50 1089 : void GREX::cmd(std::string_view key,const TypesafePtr & val) {
51 : // Enumerate all possible commands:
52 : enum {
53 : #include "GREXEnum.inc"
54 : };
55 :
56 : // Static object (initialized once) containing the map of commands:
57 : const static Tools::FastStringUnorderedMap<int> word_map = {
58 : #include "GREXMap.inc"
59 1089 : };
60 :
61 : gch::small_vector<std::string_view> words;
62 1089 : Tools::getWordsSimple(words,key);
63 1089 : unsigned nw=words.size();
64 1089 : if(nw==0) {
65 : // do nothing
66 : } else {
67 : int iword=-1;
68 : const auto it=word_map.find(words[0]);
69 1089 : if(it!=word_map.end()) {
70 1089 : iword=it->second;
71 : }
72 1089 : switch(iword) {
73 : case cmd_initialized:
74 0 : CHECK_NOTNULL(val,key);
75 0 : val.set(int(initialized));
76 : break;
77 208 : case cmd_setMPIIntracomm:
78 208 : CHECK_NOTINIT(initialized,key);
79 416 : intracomm.Set_comm(val.get<const void*>());
80 208 : break;
81 160 : case cmd_setMPIIntercomm:
82 160 : CHECK_NOTINIT(initialized,key);
83 320 : intercomm.Set_comm(val.get<const void*>());
84 320 : plumedMain.multi_sim_comm.Set_comm(val.get<const void*>());
85 160 : break;
86 0 : case cmd_setMPIFIntracomm:
87 0 : CHECK_NOTINIT(initialized,key);
88 0 : intracomm.Set_fcomm(val.get<const void*>());
89 0 : break;
90 0 : case cmd_setMPIFIntercomm:
91 0 : CHECK_NOTINIT(initialized,key);
92 0 : intercomm.Set_fcomm(val.get<const void*>());
93 0 : plumedMain.multi_sim_comm.Set_fcomm(val.get<const void*>());
94 0 : break;
95 208 : case cmd_init:
96 208 : CHECK_NOTINIT(initialized,key);
97 208 : initialized=true;
98 : // note that for PEs!=root this is automatically 0 (comm defaults to MPI_COMM_SELF)
99 208 : myreplica=intercomm.Get_rank();
100 208 : intracomm.Sum(myreplica);
101 : {
102 : std::string s;
103 208 : Tools::convert(myreplica,s);
104 416 : plumedMain.setSuffix("."+s);
105 : }
106 208 : break;
107 57 : case cmd_prepare:
108 57 : CHECK_INIT(initialized,key);
109 57 : if(intracomm.Get_rank()==0) {
110 : return;
111 : }
112 57 : intracomm.Bcast(partner,0);
113 57 : calculate();
114 : break;
115 57 : case cmd_setPartner:
116 57 : CHECK_INIT(initialized,key);
117 57 : partner=val.get<int>();
118 57 : break;
119 114 : case cmd_savePositions:
120 114 : CHECK_INIT(initialized,key);
121 114 : savePositions();
122 : break;
123 57 : case cmd_calculate:
124 57 : CHECK_INIT(initialized,key);
125 57 : if(intracomm.Get_rank()!=0) {
126 : return;
127 : }
128 57 : intracomm.Bcast(partner,0);
129 57 : calculate();
130 : break;
131 0 : case cmd_getLocalDeltaBias:
132 0 : CHECK_INIT(initialized,key);
133 0 : CHECK_NOTNULL(val,key);
134 0 : plumedMain.plumedQuantityToMD("energy",localDeltaBias,val);
135 0 : break;
136 0 : case cmd_cacheLocalUNow:
137 0 : CHECK_INIT(initialized,key);
138 0 : CHECK_NOTNULL(val,key);
139 : {
140 0 : localUNow=plumedMain.MDQuantityToPLUMED("energy",val);
141 0 : intracomm.Sum(localUNow);
142 : }
143 : break;
144 0 : case cmd_cacheLocalUSwap:
145 0 : CHECK_INIT(initialized,key);
146 0 : CHECK_NOTNULL(val,key);
147 : {
148 0 : localUSwap=plumedMain.MDQuantityToPLUMED("energy",val);
149 0 : intracomm.Sum(localUSwap);
150 : }
151 : break;
152 0 : case cmd_getForeignDeltaBias:
153 0 : CHECK_INIT(initialized,key);
154 0 : CHECK_NOTNULL(val,key);
155 0 : plumedMain.plumedQuantityToMD("energy",foreignDeltaBias,val);
156 0 : break;
157 57 : case cmd_shareAllDeltaBias:
158 57 : CHECK_INIT(initialized,key);
159 57 : if(intracomm.Get_rank()!=0) {
160 : return;
161 : }
162 57 : allDeltaBias.assign(intercomm.Get_size(),0.0);
163 57 : allDeltaBias[intercomm.Get_rank()]=localDeltaBias;
164 57 : intercomm.Sum(allDeltaBias);
165 : break;
166 171 : case cmd_getDeltaBias:
167 171 : CHECK_INIT(initialized,key);
168 171 : CHECK_NOTNULL(val,key);
169 171 : plumed_assert(nw==2);
170 171 : plumed_massert(allDeltaBias.size()==static_cast<unsigned>(intercomm.Get_size()),
171 : "to retrieve bias with cmd(\"GREX getDeltaBias\"), first share it with cmd(\"GREX shareAllDeltaBias\")");
172 : {
173 : unsigned rep;
174 171 : Tools::convert(std::string(words[1]),rep);
175 171 : plumed_massert(rep<allDeltaBias.size(),"replica index passed to cmd(\"GREX getDeltaBias\") is out of range");
176 171 : plumedMain.plumedQuantityToMD("energy",allDeltaBias[rep],val);
177 : }
178 171 : break;
179 0 : default:
180 0 : plumed_error() << "cannot interpret cmd(\" GREX" << key << "\"). check plumed developers manual to see the available commands.";
181 : break;
182 : }
183 : }
184 : }
185 :
186 114 : void GREX::savePositions() {
187 114 : plumedMain.prepareDependencies();
188 114 : plumedMain.resetActive(true);
189 114 : plumedMain.shareAll();
190 114 : plumedMain.waitData();
191 114 : std::ostringstream o;
192 114 : plumedMain.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 : plumedMain.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 : }
|