Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 2 : Copyright (c) 2016-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 "Function.h" 23 : #include "core/ActionRegister.h" 24 : #include "tools/OpenMP.h" 25 : 26 : namespace PLMD { 27 : namespace function { 28 : 29 : //+PLUMEDOC FUNCTION LOCALENSEMBLE 30 : /* 31 : Calculates the average over multiple arguments. 32 : 33 : If more than one collective variable is given for each argument then they 34 : are averaged separately. The average is stored in a component labelled _label_.cvlabel. 35 : 36 : ## Examples 37 : 38 : The following input tells plumed to calculate the chemical shifts for four 39 : different proteins in the same simulation box then average them, calculated 40 : the sum of the squared deviation with respect to the experimental values and 41 : applies a linear restraint. 42 : 43 : ```plumed 44 : MOLINFO STRUCTURE=data/template.pdb 45 : 46 : chaina: GROUP ATOMS=1-1640 47 : chainb: GROUP ATOMS=1641-3280 48 : chainc: GROUP ATOMS=3281-4920 49 : chaind: GROUP ATOMS=4921-6560 50 : 51 : WHOLEMOLECULES ENTITY0=chaina ENTITY1=chainb ENTITY2=chainc ENTITY3=chaind 52 : 53 : csa: CS2BACKBONE ATOMS=chaina NRES=100 DATA=data/ TEMPLATE=chaina.pdb NOPBC 54 : csb: CS2BACKBONE ATOMS=chainb NRES=100 DATA=data/ TEMPLATE=chainb.pdb NOPBC 55 : csc: CS2BACKBONE ATOMS=chainc NRES=100 DATA=data/ TEMPLATE=chainc.pdb NOPBC 56 : csd: CS2BACKBONE ATOMS=chaind NRES=100 DATA=data/ TEMPLATE=chaind.pdb NOPBC 57 : 58 : ensca: LOCALENSEMBLE NUM=4 ARG1=(csa\.ca_.*) ARG2=(csb\.ca_.*) ARG3=(csc\.ca_.*) ARG4=(csd\.ca_.*) 59 : enscb: LOCALENSEMBLE NUM=4 ARG1=(csa\.cb_.*) ARG2=(csb\.cb_.*) ARG3=(csc\.cb_.*) ARG4=(csd\.cb_.*) 60 : ensco: LOCALENSEMBLE NUM=4 ARG1=(csa\.co_.*) ARG2=(csb\.co_.*) ARG3=(csc\.co_.*) ARG4=(csd\.co_.*) 61 : enshn: LOCALENSEMBLE NUM=4 ARG1=(csa\.hn_.*) ARG2=(csb\.hn_.*) ARG3=(csc\.hn_.*) ARG4=(csd\.hn_.*) 62 : ensnh: LOCALENSEMBLE NUM=4 ARG1=(csa\.nh_.*) ARG2=(csb\.nh_.*) ARG3=(csc\.nh_.*) ARG4=(csd\.nh_.*) 63 : 64 : stca: STATS ARG=(ensca\.csa\.ca_.*) PARARG=(csa\.expca_.*) SQDEVSUM 65 : stcb: STATS ARG=(enscb\.csa\.cb_.*) PARARG=(csa\.expcb_.*) SQDEVSUM 66 : stco: STATS ARG=(ensco\.csa\.co_.*) PARARG=(csa\.expco_.*) SQDEVSUM 67 : sthn: STATS ARG=(enshn\.csa\.hn_.*) PARARG=(csa\.exphn_.*) SQDEVSUM 68 : stnh: STATS ARG=(ensnh\.csa\.nh_.*) PARARG=(csa\.expnh_.*) SQDEVSUM 69 : 70 : res: RESTRAINT ARG=stca.*,stcb.*,stco.*,sthn.*,stnh.* AT=0.,0.,0.,0.,0. KAPPA=0.,0.,0.,0.,0 SLOPE=16.,16.,12.,24.,0.5 71 : ``` 72 : 73 : */ 74 : //+ENDPLUMEDOC 75 : 76 : 77 : class LocalEnsemble : 78 : public Function { 79 : unsigned ens_dim; 80 : unsigned narg; 81 : public: 82 : explicit LocalEnsemble(const ActionOptions&); 83 : void calculate() override; 84 : static void registerKeywords(Keywords& keys); 85 : std::string getOutputComponentDescription( const std::string& cname, const Keywords& keys ) const override ; 86 : }; 87 : 88 : 89 : PLUMED_REGISTER_ACTION(LocalEnsemble,"LOCALENSEMBLE") 90 : 91 4 : void LocalEnsemble::registerKeywords(Keywords& keys) { 92 4 : Function::registerKeywords(keys); 93 8 : keys.addInputKeyword("numbered","ARG","scalar","the labels of the actions that you want to use"); 94 4 : keys.add("compulsory","NUM","the number of local replicas"); 95 4 : useCustomisableComponents(keys); 96 4 : } 97 : 98 2 : LocalEnsemble::LocalEnsemble(const ActionOptions&ao): 99 : Action(ao), 100 : Function(ao), 101 2 : ens_dim(0) { 102 2 : parse("NUM",ens_dim); 103 2 : if(ens_dim==0) { 104 0 : error("NUM should be greater or equal to 1"); 105 : } 106 : 107 : std::vector<Value*> arg; 108 : int oldsize=-1; 109 8 : for(unsigned i=1; i<=ens_dim; ++i ) { 110 : std::vector<Value*> larg; 111 12 : if(!parseArgumentList("ARG",i,larg)) { 112 : break; 113 : } 114 18 : for(unsigned j=0; j<larg.size(); j++) { 115 12 : arg.push_back(larg[j]); 116 : } 117 6 : if(oldsize!=-1&&oldsize!=static_cast<int>(larg.size())) { 118 0 : error("In LOCALENSEMBLE you should have the same number of arguments for each ARG keyword"); 119 : } 120 6 : oldsize = larg.size(); 121 6 : if(!larg.empty()) { 122 6 : log.printf(" with arguments %u: ", i); 123 18 : for(unsigned j=0; j<larg.size(); j++) { 124 12 : log.printf(" %s",larg[j]->getName().c_str()); 125 : } 126 6 : log.printf("\n"); 127 : } 128 : } 129 2 : requestArguments(arg); 130 2 : narg = arg.size()/ens_dim; 131 : 132 : // these are the averages 133 6 : for(unsigned i=0; i<narg; i++) { 134 4 : std::string s=getPntrToArgument(i)->getName(); 135 4 : addComponentWithDerivatives(s); 136 4 : getPntrToComponent(i)->setNotPeriodic(); 137 : } 138 : 139 2 : log.printf(" averaging over %u replicas.\n", ens_dim); 140 2 : } 141 : 142 0 : std::string LocalEnsemble::getOutputComponentDescription( const std::string& cname, const Keywords& keys ) const { 143 0 : return "the average for argument named " + cname; 144 : } 145 : 146 40 : void LocalEnsemble::calculate() { 147 40 : const double fact = 1.0/static_cast<double>(ens_dim); 148 40 : #pragma omp parallel for num_threads(OpenMP::getNumThreads()) 149 : for(unsigned i=0; i<narg; ++i) { 150 : double mean = 0.; 151 : Value* v=getPntrToComponent(i); 152 : for(unsigned j=0; j<ens_dim; ++j) { 153 : const unsigned index = j*narg+i; 154 : setDerivative(v, index, fact); 155 : mean += fact*getArgument(index); 156 : } 157 : v->set(mean); 158 : } 159 40 : } 160 : 161 : } 162 : } 163 : 164 :