GeNN  4.9.0
GPU enhanced Neuronal Networks (GeNN)
groupMerged.h
Go to the documentation of this file.
1 #pragma once
2 
3 // Standard includes
4 #include <algorithm>
5 #include <functional>
6 #include <type_traits>
7 #include <vector>
8 
9 // GeNN includes
10 #include "gennExport.h"
11 #include "currentSourceInternal.h"
12 #include "customUpdateInternal.h"
13 #include "neuronGroupInternal.h"
14 #include "synapseGroupInternal.h"
15 
16 // GeNN code generator includes
19 
20 // Forward declarations
21 namespace CodeGenerator
22 {
23 class CodeStream;
24 }
25 
26 //----------------------------------------------------------------------------
27 // CodeGenerator::GroupMerged
28 //----------------------------------------------------------------------------
30 namespace CodeGenerator
31 {
32 template<typename G>
34 {
35 public:
36  //------------------------------------------------------------------------
37  // Enumerations
38  //------------------------------------------------------------------------
39  enum class FieldType
40  {
41  Standard,
42  Host,
43  ScalarEGP,
44  PointerEGP,
45  };
46 
47  //------------------------------------------------------------------------
48  // Typedefines
49  //------------------------------------------------------------------------
50  typedef G GroupInternal;
51  typedef std::function<std::string(const G &, size_t)> GetFieldValueFunc;
52  typedef std::tuple<std::string, std::string, GetFieldValueFunc, FieldType> Field;
53 
54 
55  GroupMerged(size_t index, const std::string &precision, const std::vector<std::reference_wrapper<const GroupInternal>> groups)
56  : m_Index(index), m_LiteralSuffix((precision == "float") ? "f" : ""), m_Groups(std::move(groups))
57  {}
58 
59  //------------------------------------------------------------------------
60  // Public API
61  //------------------------------------------------------------------------
62  size_t getIndex() const { return m_Index; }
63 
65  const GroupInternal &getArchetype() const { return m_Groups.front().get(); }
66 
68  const std::string &getMemorySpace() const { return m_MemorySpace; }
69 
71  const std::vector<std::reference_wrapper<const GroupInternal>> &getGroups() const{ return m_Groups; }
72 
74  const std::vector<Field> &getFields() const{ return m_Fields; }
75 
77  std::vector<Field> getSortedFields(const BackendBase &backend) const
78  {
79  // Make a copy of fields and sort so largest come first. This should mean that due
80  // to structure packing rules, significant memory is saved and estimate is more precise
81  auto sortedFields = m_Fields;
82  std::sort(sortedFields.begin(), sortedFields.end(),
83  [&backend](const Field &a, const Field &b)
84  {
85  return (backend.getSize(std::get<0>(a)) > backend.getSize(std::get<0>(b)));
86  });
87  return sortedFields;
88 
89  }
90 
92  void generateStruct(CodeStream &os, const BackendBase &backend, const std::string &name,
93  bool host = false) const
94  {
95  os << "struct Merged" << name << "Group" << getIndex() << std::endl;
96  {
97  // Loop through fields and write to structure
98  CodeStream::Scope b(os);
99  const auto sortedFields = getSortedFields(backend);
100  for(const auto &f : sortedFields) {
101  // If field is a pointer and not marked as being a host field
102  // (in which case the backend should leave its type alone!)
103  const std::string &type = std::get<0>(f);
104  if(::Utils::isTypePointer(type) && std::get<3>(f) != FieldType::Host) {
105  // If we are generating a host structure, allow the backend to override the type
106  if(host) {
107  os << backend.getMergedGroupFieldHostType(type);
108  }
109  // Otherwise, allow the backend to add a prefix
110  else {
111  os << backend.getPointerPrefix() << type;
112  }
113  }
114  // Otherwise, leave the type alone
115  else {
116  os << type;
117  }
118  os << " " << std::get<1>(f) << ";" << std::endl;
119  }
120  os << std::endl;
121  }
122 
123  os << ";" << std::endl;
124  }
125 
127  {
128  // Get sorted fields
129  const auto sortedFields = getSortedFields(backend);
130  for(size_t fieldIndex = 0; fieldIndex < sortedFields.size(); fieldIndex++) {
131  const auto &f = sortedFields[fieldIndex];
132  os << backend.getMergedGroupFieldHostType(std::get<0>(f)) << " " << std::get<1>(f);
133  if(fieldIndex != (sortedFields.size() - 1)) {
134  os << ", ";
135  }
136  }
137  }
138 
139  size_t getStructArraySize(const BackendBase &backend) const
140  {
141  // Loop through fields again to generate any EGP pushing functions that are required and to calculate struct size
142  size_t structSize = 0;
143  size_t largestFieldSize = 0;
144  const auto sortedFields = getSortedFields(backend);
145  for(const auto &f : sortedFields) {
146  // Add size of field to total
147  const size_t fieldSize = backend.getSize(std::get<0>(f));
148  structSize += fieldSize;
149 
150  // Update largest field size
151  largestFieldSize = std::max(fieldSize, largestFieldSize);
152  }
153 
154  // Add total size of array of merged structures to merged struct data
155  // **NOTE** to match standard struct packing rules we pad to a multiple of the largest field size
156  return padSize(structSize, largestFieldSize) * getGroups().size();
157  }
158 
160 
161  void assignMemorySpaces(const BackendBase &backend, BackendBase::MemorySpaces &memorySpaces)
162  {
163  // If this backend uses memory spaces
164  if(!memorySpaces.empty()) {
165  // Get size of group in bytes
166  const size_t groupBytes = getStructArraySize(backend);
167 
168  // Loop through memory spaces
169  for(auto &m : memorySpaces) {
170  // If there is space in this memory space for group
171  if(m.second > groupBytes) {
172  // Cache memory space name in object
173  m_MemorySpace = m.first;
174 
175  // Subtract
176  m.second -= groupBytes;
177 
178  // Stop searching
179  break;
180  }
181  }
182  }
183  }
184 
185 protected:
186  //------------------------------------------------------------------------
187  // Protected methods
188  //------------------------------------------------------------------------
190  bool isParamReferenced(const std::vector<std::string> &codeStrings, const std::string &paramName) const
191  {
192  return std::any_of(codeStrings.begin(), codeStrings.end(),
193  [&paramName](const std::string &c)
194  {
195  return (c.find("$(" + paramName + ")") != std::string::npos);
196  });
197  }
198 
200  template<typename P>
201  bool isParamValueHeterogeneous(size_t index, P getParamValuesFn) const
202  {
203  // Get value of parameter in archetype group
204  const double archetypeValue = getParamValuesFn(getArchetype()).at(index);
205 
206  // Return true if any parameter values differ from the archetype value
207  return std::any_of(getGroups().cbegin(), getGroups().cend(),
208  [archetypeValue, index, getParamValuesFn](const GroupInternal &g)
209  {
210  return (getParamValuesFn(g).at(index) != archetypeValue);
211  });
212  }
213 
214  void addField(const std::string &type, const std::string &name, GetFieldValueFunc getFieldValue, FieldType fieldType = FieldType::Standard)
215  {
216  // Add field to data structure
217  m_Fields.emplace_back(type, name, getFieldValue, fieldType);
218  }
219 
220  void addScalarField(const std::string &name, GetFieldValueFunc getFieldValue, FieldType fieldType = FieldType::Standard)
221  {
222  addField("scalar", name,
223  [getFieldValue, this](const G &g, size_t i)
224  {
225  return getFieldValue(g, i) + m_LiteralSuffix;
226  },
227  fieldType);
228  }
229 
230  void addPointerField(const std::string &type, const std::string &name, const std::string &prefix)
231  {
232  assert(!Utils::isTypePointer(type));
233  addField(type + "*", name, [prefix](const G &g, size_t) { return prefix + g.getName(); });
234  }
235 
236 
237  void addVars(const Models::Base::VarVec &vars, const std::string &arrayPrefix)
238  {
239  // Loop through variables
240  for(const auto &v : vars) {
241  addPointerField(v.type, v.name, arrayPrefix + v.name);
242  }
243  }
244 
245  template<typename V>
246  void addVarReferences(const Models::Base::VarRefVec &varReferences, const std::string &arrayPrefix, V getVarRefFn)
247  {
248  // Loop through variables
249  for(size_t v = 0; v < varReferences.size(); v++) {
250  addField(varReferences[v].type + "*", varReferences[v].name,
251  [getVarRefFn, arrayPrefix, v](const G &g, size_t)
252  {
253  const auto varRef = getVarRefFn(g).at(v);
254  return arrayPrefix + varRef.getVar().name + varRef.getTargetName();
255  });
256  }
257  }
258 
259  void addEGPs(const Snippet::Base::EGPVec &egps, const std::string &arrayPrefix, const std::string &varName = "")
260  {
261  for(const auto &e : egps) {
262  const bool isPointer = Utils::isTypePointer(e.type);
263  const std::string prefix = isPointer ? arrayPrefix : "";
264  addField(e.type, e.name + varName,
265  [e, prefix, varName](const G &g, size_t) { return prefix + e.name + varName + g.getName(); },
267  }
268  }
269 
270  template<typename E>
271  void addEGPReferences(const Models::Base::EGPRefVec &egpRefs, const std::string &arrayPrefix, E getEGPRefFn)
272  {
273  for(size_t i = 0; i < egpRefs.size(); i++) {
274  const auto &e = egpRefs.at(i);
275  addField(e.type, e.name,
276  [arrayPrefix, e, getEGPRefFn, i](const G &g, size_t)
277  {
278  const auto egpRef = getEGPRefFn(g).at(i);
279  return arrayPrefix + egpRef.getEGP().name + egpRef.getTargetName();
280  },
282  }
283  }
284 
285  template<typename T, typename P, typename H>
286  void addHeterogeneousParams(const Snippet::Base::StringVec &paramNames, const std::string &suffix,
287  P getParamValues, H isHeterogeneous)
288  {
289  // Loop through params
290  for(size_t p = 0; p < paramNames.size(); p++) {
291  // If parameters is heterogeneous
292  if((static_cast<const T*>(this)->*isHeterogeneous)(p)) {
293  // Add field
294  addScalarField(paramNames[p] + suffix,
295  [p, getParamValues](const G &g, size_t)
296  {
297  const auto &values = getParamValues(g);
298  return Utils::writePreciseString(values.at(p));
299  });
300  }
301  }
302  }
303 
304  template<typename T, typename D, typename H>
305  void addHeterogeneousDerivedParams(const Snippet::Base::DerivedParamVec &derivedParams, const std::string &suffix,
306  D getDerivedParamValues, H isHeterogeneous)
307  {
308  // Loop through derived params
309  for(size_t p = 0; p < derivedParams.size(); p++) {
310  // If parameters isn't homogeneous
311  if((static_cast<const T*>(this)->*isHeterogeneous)(p)) {
312  // Add field
313  addScalarField(derivedParams[p].name + suffix,
314  [p, getDerivedParamValues](const G &g, size_t)
315  {
316  const auto &values = getDerivedParamValues(g);
317  return Utils::writePreciseString(values.at(p));
318  });
319  }
320  }
321  }
322 
323  template<typename T, typename V, typename H>
324  void addHeterogeneousVarInitParams(const Models::Base::VarVec &vars, V getVarInitialisers, H isHeterogeneous)
325  {
326  // Loop through weight update model variables
327  const std::vector<Models::VarInit> &archetypeVarInitialisers = (getArchetype().*getVarInitialisers)();
328  for(size_t v = 0; v < archetypeVarInitialisers.size(); v++) {
329  // Loop through parameters
330  const Models::VarInit &varInit = archetypeVarInitialisers[v];
331  for(size_t p = 0; p < varInit.getParams().size(); p++) {
332  if((static_cast<const T*>(this)->*isHeterogeneous)(v, p)) {
333  addScalarField(varInit.getSnippet()->getParamNames()[p] + vars[v].name,
334  [p, v, getVarInitialisers](const G &g, size_t)
335  {
336  const auto &values = (g.*getVarInitialisers)()[v].getParams();
337  return Utils::writePreciseString(values.at(p));
338  });
339  }
340  }
341  }
342  }
343 
344  template<typename T, typename V, typename H>
345  void addHeterogeneousVarInitDerivedParams(const Models::Base::VarVec &vars, V getVarInitialisers, H isHeterogeneous)
346  {
347  // Loop through weight update model variables
348  const std::vector<Models::VarInit> &archetypeVarInitialisers = (getArchetype().*getVarInitialisers)();
349  for(size_t v = 0; v < archetypeVarInitialisers.size(); v++) {
350  // Loop through parameters
351  const Models::VarInit &varInit = archetypeVarInitialisers[v];
352  for(size_t d = 0; d < varInit.getDerivedParams().size(); d++) {
353  if((static_cast<const T*>(this)->*isHeterogeneous)(v, d)) {
354  addScalarField(varInit.getSnippet()->getDerivedParams()[d].name + vars[v].name,
355  [d, v, getVarInitialisers](const G &g, size_t)
356  {
357  const auto &values = (g.*getVarInitialisers)()[v].getDerivedParams();
358  return Utils::writePreciseString(values.at(d));
359  });
360  }
361  }
362  }
363  }
364 
366  template<typename H>
367  void updateHash(H getHashableFn, boost::uuids::detail::sha1 &hash) const
368  {
369  for(const auto &g : getGroups()) {
370  Utils::updateHash(getHashableFn(g.get()), hash);
371  }
372  }
373 
374  template<typename T, typename V, typename R>
375  void updateParamHash(R isParamReferencedFn, V getValueFn, boost::uuids::detail::sha1 &hash) const
376  {
377  // Loop through parameters
378  const auto &archetypeParams = getValueFn(getArchetype());
379  for(size_t p = 0; p < archetypeParams.size(); p++) {
380  // If any of the code strings reference the parameter
381  if((static_cast<const T*>(this)->*isParamReferencedFn)(p)) {
382  // Loop through groups
383  for(const auto &g : getGroups()) {
384  // Update hash with parameter value
385  Utils::updateHash(getValueFn(g.get()).at(p), hash);
386  }
387  }
388  }
389  }
390 
391  template<typename T, typename V, typename R>
392  void updateVarInitParamHash(V getVarInitialisers, R isParamReferencedFn, boost::uuids::detail::sha1 &hash) const
393  {
394  // Loop through variables
395  const std::vector<Models::VarInit> &archetypeVarInitialisers = (getArchetype().*getVarInitialisers)();
396  for(size_t v = 0; v < archetypeVarInitialisers.size(); v++) {
397  // Loop through parameters
398  const Models::VarInit &varInit = archetypeVarInitialisers[v];
399  for(size_t p = 0; p < varInit.getParams().size(); p++) {
400  // If any of the code strings reference the parameter
401  if((static_cast<const T *>(this)->*isParamReferencedFn)(v, p)) {
402  // Loop through groups
403  for(const auto &g : getGroups()) {
404  const auto &values = (g.get().*getVarInitialisers)()[v].getParams();
405 
406  // Update hash with parameter value
407  Utils::updateHash(values.at(p), hash);
408  }
409  }
410  }
411  }
412  }
413 
414  template<typename T, typename V, typename R>
415  void updateVarInitDerivedParamHash(V getVarInitialisers, R isDerivedParamReferencedFn, boost::uuids::detail::sha1 &hash) const
416  {
417  // Loop through variables
418  const std::vector<Models::VarInit> &archetypeVarInitialisers = (getArchetype().*getVarInitialisers)();
419  for(size_t v = 0; v < archetypeVarInitialisers.size(); v++) {
420  // Loop through parameters
421  const Models::VarInit &varInit = archetypeVarInitialisers[v];
422  for(size_t d = 0; d < varInit.getDerivedParams().size(); d++) {
423  // If any of the code strings reference the parameter
424  if((static_cast<const T *>(this)->*isDerivedParamReferencedFn)(v, d)) {
425  // Loop through groups
426  for(const auto &g : getGroups()) {
427  const auto &values = (g.get().*getVarInitialisers)()[v].getDerivedParams();
428 
429  // Update hash with parameter value
430  Utils::updateHash(values.at(d), hash);
431  }
432  }
433  }
434  }
435  }
436 
437  void generateRunnerBase(const BackendBase &backend, CodeStream &definitionsInternal,
438  CodeStream &definitionsInternalFunc, CodeStream &definitionsInternalVar,
439  CodeStream &runnerVarDecl, CodeStream &runnerMergedStructAlloc,
440  const std::string &name, bool host = false) const
441  {
442  // Make a copy of fields and sort so largest come first. This should mean that due
443  // to structure packing rules, significant memory is saved and estimate is more precise
444  auto sortedFields = getSortedFields(backend);
445 
446  // If this isn't a host merged structure, generate definition for function to push group
447  if(!host) {
448  definitionsInternalFunc << "EXPORT_FUNC void pushMerged" << name << "Group" << getIndex() << "ToDevice(unsigned int idx, ";
449  generateStructFieldArgumentDefinitions(definitionsInternalFunc, backend);
450  definitionsInternalFunc << ");" << std::endl;
451  }
452 
453  // Loop through fields again to generate any EGP pushing functions that are require
454  for(const auto &f : sortedFields) {
455  // If this field is for a pointer EGP (or reference to one), also declare function to push it
456  if(std::get<3>(f) == FieldType::PointerEGP) {
457  definitionsInternalFunc << "EXPORT_FUNC void pushMerged" << name << getIndex() << std::get<1>(f) << "ToDevice(unsigned int idx, ";
458  definitionsInternalFunc << backend.getMergedGroupFieldHostType(std::get<0>(f)) << " value);" << std::endl;
459  }
460 
461  // Raise error if this field is a host field but this isn't a host structure
462  assert(std::get<3>(f) != FieldType::Host || host);
463  }
464 
465  // If merged group is used on host
466  if(host) {
467  // Generate struct directly into internal definitions
468  // **NOTE** we ignore any backend prefix as we're generating this struct for use on the host
469  generateStruct(definitionsInternal, backend, name, true);
470 
471  // Declare array of these structs containing individual neuron group pointers etc
472  runnerVarDecl << "Merged" << name << "Group" << getIndex() << " merged" << name << "Group" << getIndex() << "[" << getGroups().size() << "];" << std::endl;
473 
474  // Export it
475  definitionsInternalVar << "EXPORT_VAR Merged" << name << "Group" << getIndex() << " merged" << name << "Group" << getIndex() << "[" << getGroups().size() << "]; " << std::endl;
476  }
477 
478  // Loop through groups
479  for(size_t groupIndex = 0; groupIndex < getGroups().size(); groupIndex++) {
480  // If this is a merged group used on the host, directly set array entry
481  if(host) {
482  runnerMergedStructAlloc << "merged" << name << "Group" << getIndex() << "[" << groupIndex << "] = {";
483  generateStructFieldArguments(runnerMergedStructAlloc, groupIndex, sortedFields);
484  runnerMergedStructAlloc << "};" << std::endl;
485  }
486  // Otherwise, call function to push to device
487  else {
488  runnerMergedStructAlloc << "pushMerged" << name << "Group" << getIndex() << "ToDevice(" << groupIndex << ", ";
489  generateStructFieldArguments(runnerMergedStructAlloc, groupIndex, sortedFields);
490  runnerMergedStructAlloc << ");" << std::endl;
491  }
492  }
493  }
494 
495 private:
496  //------------------------------------------------------------------------
497  // Private methods
498  //------------------------------------------------------------------------
499  void generateStructFieldArguments(CodeStream &os, size_t groupIndex,
500  const std::vector<Field> &sortedFields) const
501  {
502  // Get group by index
503  const auto &g = getGroups()[groupIndex];
504 
505  // Loop through fields
506  for(size_t fieldIndex = 0; fieldIndex < sortedFields.size(); fieldIndex++) {
507  const auto &f = sortedFields[fieldIndex];
508  const std::string fieldInitVal = std::get<2>(f)(g, groupIndex);
509  os << fieldInitVal;
510  if(fieldIndex != (sortedFields.size() - 1)) {
511  os << ", ";
512  }
513  }
514  }
515 
516  //------------------------------------------------------------------------
517  // Members
518  //------------------------------------------------------------------------
519  const size_t m_Index;
520  const std::string m_LiteralSuffix;
521  std::string m_MemorySpace;
522  std::vector<Field> m_Fields;
523  std::vector<std::reference_wrapper<const GroupInternal>> m_Groups;
524 };
525 
526 //----------------------------------------------------------------------------
527 // CodeGenerator::NeuronSpikeQueueUpdateGroupMerged
528 //----------------------------------------------------------------------------
529 class GENN_EXPORT NeuronSpikeQueueUpdateGroupMerged : public GroupMerged<NeuronGroupInternal>
530 {
531 public:
532  NeuronSpikeQueueUpdateGroupMerged(size_t index, const std::string &precision, const std::string &timePrecison, const BackendBase &backend,
533  const std::vector<std::reference_wrapper<const NeuronGroupInternal>> &groups);
534 
535  //------------------------------------------------------------------------
536  // Public API
537  //------------------------------------------------------------------------
538  void generateRunner(const BackendBase &backend, CodeStream &definitionsInternal,
539  CodeStream &definitionsInternalFunc, CodeStream &definitionsInternalVar,
540  CodeStream &runnerVarDecl, CodeStream &runnerMergedStructAlloc) const
541  {
542  generateRunnerBase(backend, definitionsInternal, definitionsInternalFunc, definitionsInternalVar,
543  runnerVarDecl, runnerMergedStructAlloc, name);
544  }
545 
546  void genMergedGroupSpikeCountReset(CodeStream &os, unsigned int batchSize) const;
547 
548  //----------------------------------------------------------------------------
549  // Static constants
550  //----------------------------------------------------------------------------
551  static const std::string name;
552 };
553 
554 //----------------------------------------------------------------------------
555 // CodeGenerator::NeuronPrevSpikeTimeUpdateGroupMerged
556 //----------------------------------------------------------------------------
558 {
559 public:
560  NeuronPrevSpikeTimeUpdateGroupMerged(size_t index, const std::string &precision, const std::string &timePrecison, const BackendBase &backend,
561  const std::vector<std::reference_wrapper<const NeuronGroupInternal>> &groups);
562 
563  //------------------------------------------------------------------------
564  // Public API
565  //------------------------------------------------------------------------
566  void generateRunner(const BackendBase &backend, CodeStream &definitionsInternal,
567  CodeStream &definitionsInternalFunc, CodeStream &definitionsInternalVar,
568  CodeStream &runnerVarDecl, CodeStream &runnerMergedStructAlloc) const
569  {
570  generateRunnerBase(backend, definitionsInternal, definitionsInternalFunc, definitionsInternalVar,
571  runnerVarDecl, runnerMergedStructAlloc, name);
572  }
573 
574  //----------------------------------------------------------------------------
575  // Static constants
576  //----------------------------------------------------------------------------
577  static const std::string name;
578 };
579 
580 //----------------------------------------------------------------------------
581 // CodeGenerator::NeuronGroupMergedBase
582 //----------------------------------------------------------------------------
583 class GENN_EXPORT NeuronGroupMergedBase : public GroupMerged<NeuronGroupInternal>
584 {
585 public:
586  //------------------------------------------------------------------------
587  // Public API
588  //------------------------------------------------------------------------
590  bool isParamHeterogeneous(size_t index) const;
591 
593  bool isDerivedParamHeterogeneous(size_t index) const;
594 
596  bool isVarInitParamHeterogeneous(size_t varIndex, size_t paramIndex) const;
597 
599  bool isVarInitDerivedParamHeterogeneous(size_t varIndex, size_t paramIndex) const;
600 
602  bool isCurrentSourceParamHeterogeneous(size_t childIndex, size_t paramIndex) const;
603 
605  bool isCurrentSourceDerivedParamHeterogeneous(size_t childIndex, size_t paramIndex) const;
606 
608  bool isCurrentSourceVarInitParamHeterogeneous(size_t childIndex, size_t varIndex, size_t paramIndex) const;
609 
611  bool isCurrentSourceVarInitDerivedParamHeterogeneous(size_t childIndex, size_t varIndex, size_t paramIndex) const;
612 
614  bool isPSMParamHeterogeneous(size_t childIndex, size_t paramIndex) const;
615 
617  bool isPSMDerivedParamHeterogeneous(size_t childIndex, size_t varIndex) const;
618 
620  bool isPSMGlobalVarHeterogeneous(size_t childIndex, size_t paramIndex) const;
621 
623  bool isPSMVarInitParamHeterogeneous(size_t childIndex, size_t varIndex, size_t paramIndex) const;
624 
626  bool isPSMVarInitDerivedParamHeterogeneous(size_t childIndex, size_t varIndex, size_t paramIndex) const;
627 
629  const std::vector<SynapseGroupInternal*> &getSortedArchetypeMergedInSyns() const { return m_SortedMergedInSyns.front(); }
630 
632  const std::vector<SynapseGroupInternal*> &getSortedArchetypeMergedPreOutputOutSyns() const { return m_SortedMergedPreOutputOutSyns.front(); }
633 
635  const std::vector<CurrentSourceInternal*> &getSortedArchetypeCurrentSources() const { return m_SortedCurrentSources.front(); }
636 
637 protected:
638  //------------------------------------------------------------------------
639  // Protected methods
640  //------------------------------------------------------------------------
641  NeuronGroupMergedBase(size_t index, const std::string &precision, const std::string &timePrecision, const BackendBase &backend,
642  bool init, const std::vector<std::reference_wrapper<const NeuronGroupInternal>> &groups);
643 
644  void updateBaseHash(bool init, boost::uuids::detail::sha1 &hash) const;
645 
646  template<typename T, typename G, typename H>
647  void orderNeuronGroupChildren(std::vector<std::vector<T*>> &sortedGroupChildren,
648  G getVectorFunc, H getHashDigestFunc) const
649  {
650  const std::vector<T*> &archetypeChildren = (getArchetype().*getVectorFunc)();
651 
652  // Reserve vector of vectors to hold children for all neuron groups, in archetype order
653  sortedGroupChildren.reserve(getGroups().size());
654 
655  // Create temporary vector of children and their digests
656  std::vector<std::pair<boost::uuids::detail::sha1::digest_type, T*>> childDigests;
657  childDigests.reserve(archetypeChildren.size());
658 
659  // Loop through groups
660  for(const auto &g : getGroups()) {
661  // Get group children
662  const std::vector<T*> &groupChildren = (g.get().*getVectorFunc)();
663  assert(groupChildren.size() == archetypeChildren.size());
664 
665  // Loop through children and add them and their digests to vector
666  childDigests.clear();
667  for(auto *c : groupChildren) {
668  childDigests.emplace_back((c->*getHashDigestFunc)(), c);
669  }
670 
671  // Sort by digest
672  std::sort(childDigests.begin(), childDigests.end(),
673  [](const std::pair<boost::uuids::detail::sha1::digest_type, T*> &a,
674  const std::pair<boost::uuids::detail::sha1::digest_type, T*> &b)
675  {
676  return (a.first < b.first);
677  });
678 
679 
680  // Reserve vector for this group's children
681  sortedGroupChildren.emplace_back();
682  sortedGroupChildren.back().reserve(groupChildren.size());
683 
684  // Copy sorted child pointers into sortedGroupChildren
685  std::transform(childDigests.cbegin(), childDigests.cend(), std::back_inserter(sortedGroupChildren.back()),
686  [](const std::pair<boost::uuids::detail::sha1::digest_type, T*> &a){ return a.second; });
687  }
688  }
689 
691  bool isVarInitParamReferenced(size_t varIndex, size_t paramIndex) const;
692 
694  bool isVarInitDerivedParamReferenced(size_t varIndex, size_t paramIndex) const;
695 
697  bool isCurrentSourceParamReferenced(size_t childIndex, size_t paramIndex) const;
698 
700  bool isCurrentSourceDerivedParamReferenced(size_t childIndex, size_t paramIndex) const;
701 
703  bool isCurrentSourceVarInitParamReferenced(size_t childIndex, size_t varIndex, size_t paramIndex) const;
704 
706  bool isCurrentSourceVarInitDerivedParamReferenced(size_t childIndex, size_t varIndex, size_t paramIndex) const;
707 
709  bool isPSMParamReferenced(size_t childIndex, size_t paramIndex) const;
710 
712  bool isPSMDerivedParamReferenced(size_t childIndex, size_t varIndex) const;
713 
715  bool isPSMGlobalVarReferenced(size_t childIndex, size_t varIndex) const;
716 
718  bool isPSMVarInitParamReferenced(size_t childIndex, size_t varIndex, size_t paramIndex) const;
719 
721  bool isPSMVarInitDerivedParamReferenced(size_t childIndex, size_t varIndex, size_t paramIndex) const;
722 
723  template<typename T, typename G>
724  bool isChildParamValueHeterogeneous(size_t childIndex, size_t paramIndex,
725  const std::vector<std::vector<T>> &sortedGroupChildren, G getParamValuesFn) const
726  {
727  // Get value of archetype derived parameter
728  const double firstValue = getParamValuesFn(sortedGroupChildren[0][childIndex]).at(paramIndex);
729 
730  // Loop through groups within merged group
731  for(size_t i = 0; i < sortedGroupChildren.size(); i++) {
732  const auto group = sortedGroupChildren[i][childIndex];
733  if(getParamValuesFn(group).at(paramIndex) != firstValue) {
734  return true;
735  }
736  }
737 
738  return false;
739  }
740 
741  template<typename T = NeuronGroupMergedBase, typename C, typename H, typename V>
743  const std::vector<std::vector<C>> &sortedGroupChildren,
744  size_t childIndex, const std::string &prefix,
745  H isChildParamHeterogeneousFn, V getValueFn)
746  {
747  // Loop through parameters
748  for(size_t p = 0; p < paramNames.size(); p++) {
749  // If parameter is heterogeneous
750  if((static_cast<const T*>(this)->*isChildParamHeterogeneousFn)(childIndex, p)) {
751  addScalarField(paramNames[p] + prefix + std::to_string(childIndex),
752  [&sortedGroupChildren, childIndex, p, getValueFn](const NeuronGroupInternal &, size_t groupIndex)
753  {
754  const auto *child = sortedGroupChildren.at(groupIndex).at(childIndex);
755  return Utils::writePreciseString((child->*getValueFn)().at(p));
756  });
757  }
758  }
759  }
760 
761  template<typename T = NeuronGroupMergedBase, typename C, typename H, typename V>
763  const std::vector<std::vector<C>> &sortedGroupChildren,
764  size_t childIndex, const std::string &prefix,
765  H isChildDerivedParamHeterogeneousFn, V getValueFn)
766  {
767  // Loop through derived parameters
768  for(size_t p = 0; p < derivedParams.size(); p++) {
769  // If parameter is heterogeneous
770  if((static_cast<const T*>(this)->*isChildDerivedParamHeterogeneousFn)(childIndex, p)) {
771  addScalarField(derivedParams[p].name + prefix + std::to_string(childIndex),
772  [&sortedGroupChildren, childIndex, p, getValueFn](const NeuronGroupInternal &, size_t groupIndex)
773  {
774  const auto *child = sortedGroupChildren.at(groupIndex).at(childIndex);
775  return Utils::writePreciseString((child->*getValueFn)().at(p));
776  });
777  }
778  }
779  }
780 
781  template<typename T = NeuronGroupMergedBase, typename C, typename H, typename V>
783  const std::vector<std::vector<C>> &sortedGroupChildren,
784  size_t childIndex, size_t varIndex, const std::string &prefix,
785  H isChildParamHeterogeneousFn, V getVarInitialiserFn)
786  {
787  // Loop through parameters
788  for(size_t p = 0; p < paramNames.size(); p++) {
789  // If parameter is heterogeneous
790  if((static_cast<const T*>(this)->*isChildParamHeterogeneousFn)(childIndex, varIndex, p)) {
791  addScalarField(paramNames[p] + prefix + std::to_string(childIndex),
792  [&sortedGroupChildren, childIndex, varIndex, p, getVarInitialiserFn](const NeuronGroupInternal &, size_t groupIndex)
793  {
794  const auto *child = sortedGroupChildren.at(groupIndex).at(childIndex);
795  const std::vector<Models::VarInit> &varInit = (child->*getVarInitialiserFn)();
796  return Utils::writePreciseString(varInit.at(varIndex).getParams().at(p));
797  });
798  }
799  }
800  }
801 
802  template<typename T = NeuronGroupMergedBase, typename C, typename H, typename V>
804  const std::vector<std::vector<C>> &sortedGroupChildren,
805  size_t childIndex, size_t varIndex, const std::string &prefix,
806  H isChildDerivedParamHeterogeneousFn, V getVarInitialiserFn)
807  {
808  // Loop through parameters
809  for(size_t p = 0; p < derivedParams.size(); p++) {
810  // If parameter is heterogeneous
811  if((static_cast<const T*>(this)->*isChildDerivedParamHeterogeneousFn)(childIndex, varIndex, p)) {
812  addScalarField(derivedParams[p].name + prefix + std::to_string(childIndex),
813  [&sortedGroupChildren, childIndex, varIndex, p, getVarInitialiserFn](const NeuronGroupInternal &, size_t groupIndex)
814  {
815  const auto *child = sortedGroupChildren.at(groupIndex).at(childIndex);
816  const std::vector<Models::VarInit> &varInit = (child->*getVarInitialiserFn)();
817  return Utils::writePreciseString(varInit.at(varIndex).getDerivedParams().at(p));
818  });
819  }
820  }
821  }
822 
823  template<typename S>
824  void addChildEGPs(const std::vector<Snippet::Base::EGP> &egps, size_t childIndex,
825  const std::string &arrayPrefix, const std::string &prefix,
826  S getEGPSuffixFn)
827  {
828  for(const auto &e : egps) {
829  const bool isPointer = Utils::isTypePointer(e.type);
830  const std::string varPrefix = isPointer ? arrayPrefix : "";
831  addField(e.type, e.name + prefix + std::to_string(childIndex),
832  [getEGPSuffixFn, childIndex, e, varPrefix](const NeuronGroupInternal&, size_t groupIndex)
833  {
834  return varPrefix + e.name + getEGPSuffixFn(groupIndex, childIndex);
835  },
837  }
838  }
839 
840  template<typename T = NeuronGroupMergedBase, typename C, typename V, typename R>
841  void updateChildParamHash(const std::vector<std::vector<C>> &sortedGroupChildren,
842  size_t childIndex, R isChildParamReferencedFn, V getValueFn,
843  boost::uuids::detail::sha1 &hash) const
844  {
845  // Loop through parameters
846  const auto &archetypeParamNames = (sortedGroupChildren.front().at(childIndex)->*getValueFn)();
847  for(size_t p = 0; p < archetypeParamNames.size(); p++) {
848  // If any of the code strings reference the parameter
849  if((static_cast<const T*>(this)->*isChildParamReferencedFn)(childIndex, p)) {
850  // Loop through groups
851  for(size_t g = 0; g < getGroups().size(); g++) {
852  // Get child group
853  const auto *child = sortedGroupChildren.at(g).at(childIndex);
854 
855  // Update hash with parameter value
856  Utils::updateHash((child->*getValueFn)().at(p), hash);
857  }
858  }
859  }
860  }
861 
862  template<typename T = NeuronGroupMergedBase, typename C, typename V, typename R>
863  void updateChildDerivedParamHash(const std::vector<std::vector<C>> &sortedGroupChildren,
864  size_t childIndex, R isChildDerivedParamReferencedFn, V getValueFn,
865  boost::uuids::detail::sha1 &hash) const
866  {
867  // Loop through derived parameters
868  const auto &archetypeDerivedParams = (sortedGroupChildren.front().at(childIndex)->*getValueFn)();
869  for(size_t p = 0; p < archetypeDerivedParams.size(); p++) {
870  // If any of the code strings reference the parameter
871  if((static_cast<const T*>(this)->*isChildDerivedParamReferencedFn)(childIndex, p)) {
872  // Loop through groups
873  for(size_t g = 0; g < getGroups().size(); g++) {
874  // Get child group
875  const auto *child = sortedGroupChildren.at(g).at(childIndex);
876 
877  // Update hash with parameter value
878  Utils::updateHash((child->*getValueFn)().at(p), hash);
879  }
880  }
881  }
882  }
883 
884  template<typename T = NeuronGroupMergedBase, typename C, typename R, typename V>
885  void updateChildVarInitParamsHash(const std::vector<std::vector<C>> &sortedGroupChildren,
886  size_t childIndex, size_t varIndex, R isChildParamReferencedFn, V getVarInitialiserFn,
887  boost::uuids::detail::sha1 &hash) const
888  {
889  // Loop through parameters
890  const auto &archetypeVarInit = (sortedGroupChildren.front().at(childIndex)->*getVarInitialiserFn)();
891  const auto &archetypeParams = archetypeVarInit.at(varIndex).getParams();
892  for(size_t p = 0; p < archetypeParams.size(); p++) {
893  // If parameter is referenced
894  if((static_cast<const T*>(this)->*isChildParamReferencedFn)(childIndex, varIndex, p)) {
895  // Loop through groups
896  for(size_t g = 0; g < getGroups().size(); g++) {
897  // Get child group and its variable initialisers
898  const auto *child = sortedGroupChildren.at(g).at(childIndex);
899  const std::vector<Models::VarInit> &varInit = (child->*getVarInitialiserFn)();
900 
901  // Update hash with parameter value
902  Utils::updateHash(varInit.at(varIndex).getParams().at(p), hash);
903  }
904  }
905  }
906  }
907 
908  template<typename T = NeuronGroupMergedBase, typename C, typename R, typename V>
909  void updateChildVarInitDerivedParamsHash(const std::vector<std::vector<C>> &sortedGroupChildren,
910  size_t childIndex, size_t varIndex, R isChildDerivedParamReferencedFn, V getVarInitialiserFn,
911  boost::uuids::detail::sha1 &hash) const
912  {
913  // Loop through derived parameters
914  const auto &archetypeVarInit = (sortedGroupChildren.front().at(childIndex)->*getVarInitialiserFn)();
915  const auto &archetypeDerivedParams = archetypeVarInit.at(varIndex).getDerivedParams();
916  for(size_t p = 0; p < archetypeDerivedParams.size(); p++) {
917  // If parameter is referenced
918  if((static_cast<const T*>(this)->*isChildDerivedParamReferencedFn)(childIndex, varIndex, p)) {
919  // Loop through groups
920  for(size_t g = 0; g < getGroups().size(); g++) {
921  // Get child group and its variable initialisers
922  const auto *child = sortedGroupChildren.at(g).at(childIndex);
923  const auto &varInit = (child->*getVarInitialiserFn)();
924 
925  // Update hash with parameter value
926  Utils::updateHash(varInit.at(varIndex).getDerivedParams().at(p), hash);
927  }
928  }
929  }
930  }
931 
932  void addMergedInSynPointerField(const std::string &type, const std::string &name,
933  size_t archetypeIndex, const std::string &prefix);
934 
935  void addMergedPreOutputOutSynPointerField(const std::string &type, const std::string &name,
936  size_t archetypeIndex, const std::string &prefix);
937 
938 
939 private:
940  //------------------------------------------------------------------------
941  // Members
942  //------------------------------------------------------------------------
943  std::vector<std::vector<SynapseGroupInternal*>> m_SortedMergedInSyns;
944  std::vector<std::vector<SynapseGroupInternal*>> m_SortedMergedPreOutputOutSyns;
945  std::vector<std::vector<CurrentSourceInternal*>> m_SortedCurrentSources;
946 };
947 
948 
949 
950 //----------------------------------------------------------------------------
951 // CodeGenerator::SynapseDendriticDelayUpdateGroupMerged
952 //----------------------------------------------------------------------------
954 {
955 public:
956  SynapseDendriticDelayUpdateGroupMerged(size_t index, const std::string &precision, const std::string &timePrecision, const BackendBase &backend,
957  const std::vector<std::reference_wrapper<const SynapseGroupInternal>> &group);
958 
959  //------------------------------------------------------------------------
960  // Public API
961  //------------------------------------------------------------------------
962  void generateRunner(const BackendBase &backend, CodeStream &definitionsInternal,
963  CodeStream &definitionsInternalFunc, CodeStream &definitionsInternalVar,
964  CodeStream &runnerVarDecl, CodeStream &runnerMergedStructAlloc) const
965  {
966  generateRunnerBase(backend, definitionsInternal, definitionsInternalFunc, definitionsInternalVar,
967  runnerVarDecl, runnerMergedStructAlloc, name);
968  }
969 
970  //----------------------------------------------------------------------------
971  // Static constants
972  //----------------------------------------------------------------------------
973  static const std::string name;
974 };
975 
976 // ----------------------------------------------------------------------------
977 // SynapseConnectivityHostInitGroupMerged
978 //----------------------------------------------------------------------------
980 {
981 public:
982  SynapseConnectivityHostInitGroupMerged(size_t index, const std::string &precision, const std::string &timePrecision, const BackendBase &backend,
983  const std::vector<std::reference_wrapper<const SynapseGroupInternal>> &groups);
984 
985  //------------------------------------------------------------------------
986  // Public API
987  //------------------------------------------------------------------------
988  void generateRunner(const BackendBase &backend, CodeStream &definitionsInternal,
989  CodeStream &definitionsInternalFunc, CodeStream &definitionsInternalVar,
990  CodeStream &runnerVarDecl, CodeStream &runnerMergedStructAlloc) const
991  {
992  generateRunnerBase(backend, definitionsInternal, definitionsInternalFunc, definitionsInternalVar,
993  runnerVarDecl, runnerMergedStructAlloc, name, true);
994  }
995 
997  bool isConnectivityInitParamHeterogeneous(size_t paramIndex) const;
998 
1000  bool isConnectivityInitDerivedParamHeterogeneous(size_t paramIndex) const;
1001 
1002  //----------------------------------------------------------------------------
1003  // Static constants
1004  //----------------------------------------------------------------------------
1005  static const std::string name;
1006 
1007 private:
1008  //------------------------------------------------------------------------
1009  // Private methods
1010  //------------------------------------------------------------------------
1012  bool isSparseConnectivityInitParamReferenced(size_t paramIndex) const;
1013 
1015  bool isSparseConnectivityInitDerivedParamReferenced(size_t paramIndex) const;
1016 };
1017 
1018 //----------------------------------------------------------------------------
1019 // CodeGenerator::SynapseGroupMergedBase
1020 //----------------------------------------------------------------------------
1021 class GENN_EXPORT SynapseGroupMergedBase : public GroupMerged<SynapseGroupInternal>
1022 {
1023 public:
1024  //------------------------------------------------------------------------
1025  // Public API
1026  //------------------------------------------------------------------------
1028  bool isWUParamHeterogeneous(size_t paramIndex) const;
1029 
1031  bool isWUDerivedParamHeterogeneous(size_t paramIndex) const;
1032 
1034  bool isWUGlobalVarHeterogeneous(size_t varIndex) const;
1035 
1037  bool isWUVarInitParamHeterogeneous(size_t varIndex, size_t paramIndex) const;
1038 
1040  bool isWUVarInitDerivedParamHeterogeneous(size_t varIndex, size_t paramIndex) const;
1041 
1043  bool isSparseConnectivityInitParamHeterogeneous(size_t paramIndex) const;
1044 
1046  bool isSparseConnectivityInitDerivedParamHeterogeneous(size_t paramIndex) const;
1047 
1049  bool isToeplitzConnectivityInitParamHeterogeneous(size_t paramIndex) const;
1050 
1052  bool isToeplitzConnectivityInitDerivedParamHeterogeneous(size_t paramIndex) const;
1053 
1055  bool isSrcNeuronParamHeterogeneous(size_t paramIndex) const;
1056 
1058  bool isSrcNeuronDerivedParamHeterogeneous(size_t paramIndex) const;
1059 
1061  bool isTrgNeuronParamHeterogeneous(size_t paramIndex) const;
1062 
1064  bool isTrgNeuronDerivedParamHeterogeneous(size_t paramIndex) const;
1065 
1067  bool isKernelSizeHeterogeneous(size_t dimensionIndex) const
1068  {
1069  return CodeGenerator::isKernelSizeHeterogeneous(this, dimensionIndex, getGroupKernelSize);
1070  }
1071 
1073  std::string getKernelSize(size_t dimensionIndex) const
1074  {
1075  return CodeGenerator::getKernelSize(this, dimensionIndex, getGroupKernelSize);
1076  }
1077 
1079  void genKernelIndex(std::ostream& os, const CodeGenerator::Substitutions& subs) const
1080  {
1081  return CodeGenerator::genKernelIndex(this, os, subs, getGroupKernelSize);
1082  }
1083 
1084  std::string getPreSlot(unsigned int batchSize) const;
1085  std::string getPostSlot(unsigned int batchSize) const;
1086 
1087  std::string getPreVarIndex(unsigned int batchSize, VarAccessDuplication varDuplication, const std::string &index) const
1088  {
1089  return getPreVarIndex(getArchetype().getSrcNeuronGroup()->isDelayRequired(), batchSize, varDuplication, index);
1090  }
1091 
1092  std::string getPostVarIndex(unsigned int batchSize, VarAccessDuplication varDuplication, const std::string &index) const
1093  {
1094  return getPostVarIndex(getArchetype().getTrgNeuronGroup()->isDelayRequired(), batchSize, varDuplication, index);
1095  }
1096 
1097  std::string getPreWUVarIndex(unsigned int batchSize, VarAccessDuplication varDuplication, const std::string &index) const
1098  {
1099  return getPreVarIndex(getArchetype().getDelaySteps() != 0, batchSize, varDuplication, index);
1100  }
1101 
1102  std::string getPostWUVarIndex(unsigned int batchSize, VarAccessDuplication varDuplication, const std::string &index) const
1103  {
1104  return getPostVarIndex(getArchetype().getBackPropDelaySteps() != 0, batchSize, varDuplication, index);
1105  }
1106 
1107  std::string getPostDenDelayIndex(unsigned int batchSize, const std::string &index, const std::string &offset) const;
1108 
1109  std::string getPreVarIndex(bool delay, unsigned int batchSize, VarAccessDuplication varDuplication, const std::string &index) const;
1110  std::string getPostVarIndex(bool delay, unsigned int batchSize, VarAccessDuplication varDuplication, const std::string &index) const;
1111 
1112  std::string getPrePrevSpikeTimeIndex(bool delay, unsigned int batchSize, VarAccessDuplication varDuplication, const std::string &index) const;
1113  std::string getPostPrevSpikeTimeIndex(bool delay, unsigned int batchSize, VarAccessDuplication varDuplication, const std::string &index) const;
1114 
1115  std::string getPostISynIndex(unsigned int batchSize, const std::string &index) const
1116  {
1117  return ((batchSize == 1) ? "" : "postBatchOffset + ") + index;
1118  }
1119 
1120  std::string getPreISynIndex(unsigned int batchSize, const std::string &index) const
1121  {
1122  return ((batchSize == 1) ? "" : "preBatchOffset + ") + index;
1123  }
1124 
1125  std::string getSynVarIndex(unsigned int batchSize, VarAccessDuplication varDuplication, const std::string &index) const;
1126  std::string getKernelVarIndex(unsigned int batchSize, VarAccessDuplication varDuplication, const std::string &index) const;
1127 
1128 protected:
1129  //----------------------------------------------------------------------------
1130  // Enumerations
1131  //----------------------------------------------------------------------------
1132  enum class Role
1133  {
1134  PresynapticUpdate,
1135  PostsynapticUpdate,
1136  SynapseDynamics,
1137  Init,
1138  SparseInit,
1139  ConnectivityInit,
1140  };
1141 
1142  SynapseGroupMergedBase(size_t index, const std::string &precision, const std::string &timePrecision, const BackendBase &backend,
1143  Role role, const std::string &archetypeCode, const std::vector<std::reference_wrapper<const SynapseGroupInternal>> &groups);
1144 
1145  //----------------------------------------------------------------------------
1146  // Protected methods
1147  //----------------------------------------------------------------------------
1148  boost::uuids::detail::sha1::digest_type getHashDigest(Role role) const;
1149 
1150  const std::string &getArchetypeCode() const { return m_ArchetypeCode; }
1151 
1152 private:
1153  //------------------------------------------------------------------------
1154  // Private methods
1155  //------------------------------------------------------------------------
1156  void addPSPointerField(const std::string &type, const std::string &name, const std::string &prefix);
1157  void addPreOutputPointerField(const std::string &type, const std::string &name, const std::string &prefix);
1158  void addSrcPointerField(const std::string &type, const std::string &name, const std::string &prefix);
1159  void addTrgPointerField(const std::string &type, const std::string &name, const std::string &prefix);
1160  void addWeightSharingPointerField(const std::string &type, const std::string &name, const std::string &prefix);
1161 
1162  std::string getVarIndex(bool delay, unsigned int batchSize, VarAccessDuplication varDuplication,
1163  const std::string &index, const std::string &prefix) const;
1164 
1166  bool isWUParamReferenced(size_t paramIndex) const;
1167 
1169  bool isWUDerivedParamReferenced(size_t paramIndex) const;
1170 
1172  bool isWUGlobalVarReferenced(size_t varIndex) const;
1173 
1175  bool isWUVarInitParamReferenced(size_t varIndex, size_t paramIndex) const;
1176 
1178  bool isWUVarInitDerivedParamReferenced(size_t varIndex, size_t paramIndex) const;
1179 
1181  bool isSparseConnectivityInitParamReferenced(size_t paramIndex) const;
1182 
1184  bool isSparseConnectivityInitDerivedParamReferenced(size_t paramIndex) const;
1185 
1187  bool isToeplitzConnectivityInitParamReferenced(size_t paramIndex) const;
1188 
1190  bool isToeplitzConnectivityInitDerivedParamReferenced(size_t paramIndex) const;
1191 
1193  bool isSrcNeuronParamReferenced(size_t paramIndex) const;
1194 
1196  bool isSrcNeuronDerivedParamReferenced(size_t paramIndex) const;
1197 
1199  bool isTrgNeuronParamReferenced(size_t paramIndex) const;
1200 
1202  bool isTrgNeuronDerivedParamReferenced(size_t paramIndex) const;
1203 
1204  static const std::vector<unsigned int>& getGroupKernelSize(const SynapseGroupInternal& g)
1205  {
1206  return g.getKernelSize();
1207  }
1208 
1209  //------------------------------------------------------------------------
1210  // Members
1211  //------------------------------------------------------------------------
1212  const std::string m_ArchetypeCode;
1213 };
1214 
1215 // ----------------------------------------------------------------------------
1216 // CustomUpdateHostReductionGroupMergedBase
1217 //----------------------------------------------------------------------------
1218 template<typename G>
1220 {
1221 protected:
1222  CustomUpdateHostReductionGroupMergedBase(size_t index, const std::string &precision, const BackendBase &backend,
1223  const std::vector<std::reference_wrapper<const G>> &groups)
1224  : GroupMerged<G>(index, precision, groups)
1225  {
1226  // Loop through variables and add pointers if they are reduction targets
1227  const CustomUpdateModels::Base *cm = this->getArchetype().getCustomUpdateModel();
1228  for(const auto &v : cm->getVars()) {
1229  if(v.access & VarAccessModeAttribute::REDUCE) {
1230  this->addPointerField(v.type, v.name, backend.getDeviceVarPrefix() + v.name);
1231  }
1232  }
1233 
1234  // Loop through variable references and add pointers if they are reduction targets
1235  for(const auto &v : cm->getVarRefs()) {
1236  if(v.access & VarAccessModeAttribute::REDUCE) {
1237  this->addPointerField(v.type, v.name, backend.getDeviceVarPrefix() + v.name);
1238  }
1239  }
1240  }
1241 };
1242 
1243 // ----------------------------------------------------------------------------
1244 // CustomUpdateHostReductionGroupMerged
1245 //----------------------------------------------------------------------------
1247 {
1248 public:
1249  CustomUpdateHostReductionGroupMerged(size_t index, const std::string &precision, const std::string &, const BackendBase &backend,
1250  const std::vector<std::reference_wrapper<const CustomUpdateInternal>> &groups);
1251 
1252  //------------------------------------------------------------------------
1253  // Public API
1254  //------------------------------------------------------------------------
1255  void generateRunner(const BackendBase &backend, CodeStream &definitionsInternal,
1256  CodeStream &definitionsInternalFunc, CodeStream &definitionsInternalVar,
1257  CodeStream &runnerVarDecl, CodeStream &runnerMergedStructAlloc) const
1258  {
1259  generateRunnerBase(backend, definitionsInternal, definitionsInternalFunc, definitionsInternalVar,
1260  runnerVarDecl, runnerMergedStructAlloc, name, true);
1261  }
1262 
1263  //----------------------------------------------------------------------------
1264  // Static constants
1265  //----------------------------------------------------------------------------
1266  static const std::string name;
1267 };
1268 
1269 // ----------------------------------------------------------------------------
1270 // CustomWUUpdateHostReductionGroupMerged
1271 //----------------------------------------------------------------------------
1273 {
1274 public:
1275  CustomWUUpdateHostReductionGroupMerged(size_t index, const std::string &precision, const std::string &, const BackendBase &backend,
1276  const std::vector<std::reference_wrapper<const CustomUpdateWUInternal>> &groups);
1277 
1278  //------------------------------------------------------------------------
1279  // Public API
1280  //------------------------------------------------------------------------
1281  void generateRunner(const BackendBase &backend, CodeStream &definitionsInternal,
1282  CodeStream &definitionsInternalFunc, CodeStream &definitionsInternalVar,
1283  CodeStream &runnerVarDecl, CodeStream &runnerMergedStructAlloc) const
1284  {
1285  generateRunnerBase(backend, definitionsInternal, definitionsInternalFunc, definitionsInternalVar,
1286  runnerVarDecl, runnerMergedStructAlloc, name, true);
1287  }
1288 
1289  //----------------------------------------------------------------------------
1290  // Static constants
1291  //----------------------------------------------------------------------------
1292  static const std::string name;
1293 };
1294 } // namespace CodeGenerator
const std::vector< SynapseGroupInternal * > & getSortedArchetypeMergedPreOutputOutSyns() const
Get sorted vectors of merged outgoing synapse groups with presynaptic output belonging to archetype g...
Definition: groupMerged.h:632
void updateChildParamHash(const std::vector< std::vector< C >> &sortedGroupChildren, size_t childIndex, R isChildParamReferencedFn, V getValueFn, boost::uuids::detail::sha1 &hash) const
Definition: groupMerged.h:841
Definition: neuronGroupInternal.h:9
void genKernelIndex(std::ostream &os, const CodeGenerator::Substitutions &subs) const
Generate an index into a kernel based on the id_kernel_XXX variables in subs.
Definition: groupMerged.h:1079
static const std::string name
Definition: groupMerged.h:973
static const std::string name
Definition: groupMerged.h:1266
std::vector< EGPRef > EGPRefVec
Definition: models.h:117
void generateRunner(const BackendBase &backend, CodeStream &definitionsInternal, CodeStream &definitionsInternalFunc, CodeStream &definitionsInternalVar, CodeStream &runnerVarDecl, CodeStream &runnerMergedStructAlloc) const
Definition: groupMerged.h:1255
std::vector< Var > VarVec
Definition: models.h:115
Definition: groupMerged.h:1021
const std::vector< CurrentSourceInternal * > & getSortedArchetypeCurrentSources() const
Get sorted vectors of current sources belonging to archetype group.
Definition: groupMerged.h:635
void addHeterogeneousChildDerivedParams(const Snippet::Base::DerivedParamVec &derivedParams, const std::vector< std::vector< C >> &sortedGroupChildren, size_t childIndex, const std::string &prefix, H isChildDerivedParamHeterogeneousFn, V getValueFn)
Definition: groupMerged.h:762
const std::vector< Field > & getFields() const
Get group fields.
Definition: groupMerged.h:74
void generateRunner(const BackendBase &backend, CodeStream &definitionsInternal, CodeStream &definitionsInternalFunc, CodeStream &definitionsInternalVar, CodeStream &runnerVarDecl, CodeStream &runnerMergedStructAlloc) const
Definition: groupMerged.h:538
void generateRunner(const BackendBase &backend, CodeStream &definitionsInternal, CodeStream &definitionsInternalFunc, CodeStream &definitionsInternalVar, CodeStream &runnerVarDecl, CodeStream &runnerMergedStructAlloc) const
Definition: groupMerged.h:962
void genKernelIndex(const G *group, std::ostream &os, const CodeGenerator::Substitutions &subs, K getKernelSizeFn)
Definition: codeGenUtils.h:180
void addHeterogeneousChildParams(const Snippet::Base::StringVec &paramNames, const std::vector< std::vector< C >> &sortedGroupChildren, size_t childIndex, const std::string &prefix, H isChildParamHeterogeneousFn, V getValueFn)
Definition: groupMerged.h:742
#define GENN_EXPORT
Definition: gennExport.h:13
void addHeterogeneousChildVarInitDerivedParams(const Snippet::Base::DerivedParamVec &derivedParams, const std::vector< std::vector< C >> &sortedGroupChildren, size_t childIndex, size_t varIndex, const std::string &prefix, H isChildDerivedParamHeterogeneousFn, V getVarInitialiserFn)
Definition: groupMerged.h:803
void addPointerField(const std::string &type, const std::string &name, const std::string &prefix)
Definition: groupMerged.h:230
void addHeterogeneousVarInitParams(const Models::Base::VarVec &vars, V getVarInitialisers, H isHeterogeneous)
Definition: groupMerged.h:324
STL namespace.
void addHeterogeneousVarInitDerivedParams(const Models::Base::VarVec &vars, V getVarInitialisers, H isHeterogeneous)
Definition: groupMerged.h:345
void generateStruct(CodeStream &os, const BackendBase &backend, const std::string &name, bool host=false) const
Generate declaration of struct to hold this merged group.
Definition: groupMerged.h:92
void addChildEGPs(const std::vector< Snippet::Base::EGP > &egps, size_t childIndex, const std::string &arrayPrefix, const std::string &prefix, S getEGPSuffixFn)
Definition: groupMerged.h:824
static const std::string name
Definition: groupMerged.h:1005
GENN_EXPORT void init(plog::Severity gennLevel, plog::Severity codeGeneratorLevel, plog::IAppender *gennAppender, plog::IAppender *codeGeneratorAppender)
Definition: logging.cc:6
void addField(const std::string &type, const std::string &name, GetFieldValueFunc getFieldValue, FieldType fieldType=FieldType::Standard)
Definition: groupMerged.h:214
void writePreciseString(std::ostream &os, T value)
This function writes a floating point value to a stream -setting the precision so no digits are lost...
Definition: gennUtils.h:92
Helper class for generating code - automatically inserts brackets, indents etc.
Definition: backendBase.h:30
void generateRunner(const BackendBase &backend, CodeStream &definitionsInternal, CodeStream &definitionsInternalFunc, CodeStream &definitionsInternalVar, CodeStream &runnerVarDecl, CodeStream &runnerMergedStructAlloc) const
Definition: groupMerged.h:566
std::string getPreISynIndex(unsigned int batchSize, const std::string &index) const
Definition: groupMerged.h:1120
void assignMemorySpaces(const BackendBase &backend, BackendBase::MemorySpaces &memorySpaces)
Assign memory spaces to group.
Definition: groupMerged.h:161
Base class for all current source models.
Definition: customUpdateModels.h:31
bool isParamReferenced(const std::vector< std::string > &codeStrings, const std::string &paramName) const
Helper to test whether parameter is referenced in vector of codestrings.
Definition: groupMerged.h:190
void updateChildVarInitParamsHash(const std::vector< std::vector< C >> &sortedGroupChildren, size_t childIndex, size_t varIndex, R isChildParamReferencedFn, V getVarInitialiserFn, boost::uuids::detail::sha1 &hash) const
Definition: groupMerged.h:885
const std::vector< std::reference_wrapper< const GroupInternal > > & getGroups() const
Gets access to underlying vector of neuron groups which have been merged.
Definition: groupMerged.h:71
Definition: groupMerged.h:33
std::vector< std::pair< std::string, size_t > > MemorySpaces
Vector of prefixes required to allocate in memory space and size of memory space. ...
Definition: backendBase.h:190
Definition: synapseGroupInternal.h:9
Definition: codeStream.h:21
const std::string & getArchetypeCode() const
Definition: groupMerged.h:1150
std::function< std::string(const G &, size_t)> GetFieldValueFunc
Definition: groupMerged.h:51
Definition: substitutions.h:21
void addHeterogeneousDerivedParams(const Snippet::Base::DerivedParamVec &derivedParams, const std::string &suffix, D getDerivedParamValues, H isHeterogeneous)
Definition: groupMerged.h:305
const std::string & getMemorySpace() const
Get name of memory space assigned to group.
Definition: groupMerged.h:68
void updateVarInitParamHash(V getVarInitialisers, R isParamReferencedFn, boost::uuids::detail::sha1 &hash) const
Definition: groupMerged.h:392
void addHeterogeneousParams(const Snippet::Base::StringVec &paramNames, const std::string &suffix, P getParamValues, H isHeterogeneous)
Definition: groupMerged.h:286
const std::vector< unsigned int > & getKernelSize() const
Definition: synapseGroup.h:130
void updateHash(const T &value, boost::uuids::detail::sha1 &hash)
Hash arithmetic types and enums.
Definition: gennUtils.h:128
size_t getIndex() const
Definition: groupMerged.h:62
bool isChildParamValueHeterogeneous(size_t childIndex, size_t paramIndex, const std::vector< std::vector< T >> &sortedGroupChildren, G getParamValuesFn) const
Definition: groupMerged.h:724
Definition: backendBase.h:176
size_t padSize(size_t size, size_t blockSize)
Pad an integer to a multiple of another.
Definition: codeGenUtils.h:64
void generateRunner(const BackendBase &backend, CodeStream &definitionsInternal, CodeStream &definitionsInternalFunc, CodeStream &definitionsInternalVar, CodeStream &runnerVarDecl, CodeStream &runnerMergedStructAlloc) const
Definition: groupMerged.h:1281
GENN_EXPORT bool isTypePointer(const std::string &type)
Function to determine whether a string containing a type is a pointer.
Definition: gennUtils.cc:75
std::vector< Field > getSortedFields(const BackendBase &backend) const
Get group fields, sorted into order they will appear in struct.
Definition: groupMerged.h:77
Definition: models.h:151
void addEGPReferences(const Models::Base::EGPRefVec &egpRefs, const std::string &arrayPrefix, E getEGPRefFn)
Definition: groupMerged.h:271
bool isKernelSizeHeterogeneous(const G *group, size_t dimensionIndex, K getKernelSizeFn)
Definition: codeGenUtils.h:152
void updateChildDerivedParamHash(const std::vector< std::vector< C >> &sortedGroupChildren, size_t childIndex, R isChildDerivedParamReferencedFn, V getValueFn, boost::uuids::detail::sha1 &hash) const
Definition: groupMerged.h:863
Definition: groupMerged.h:583
virtual StringVec getParamNames() const
Gets names of of (independent) model parameters.
Definition: snippet.h:187
const GroupInternal & getArchetype() const
Get &#39;archetype&#39; neuron group - it&#39;s properties represent those of all other merged neuron groups...
Definition: groupMerged.h:65
static const std::string name
Definition: groupMerged.h:551
static const std::string name
Definition: groupMerged.h:577
virtual VarVec getVars() const
Gets names and types (as strings) of model variables.
Definition: models.h:123
const SnippetBase * getSnippet() const
Definition: snippet.h:263
void updateVarInitDerivedParamHash(V getVarInitialisers, R isDerivedParamReferencedFn, boost::uuids::detail::sha1 &hash) const
Definition: groupMerged.h:415
GroupMerged(size_t index, const std::string &precision, const std::vector< std::reference_wrapper< const GroupInternal >> groups)
Definition: groupMerged.h:55
const std::vector< double > & getParams() const
Definition: snippet.h:264
void generateRunner(const BackendBase &backend, CodeStream &definitionsInternal, CodeStream &definitionsInternalFunc, CodeStream &definitionsInternalVar, CodeStream &runnerVarDecl, CodeStream &runnerMergedStructAlloc) const
Definition: groupMerged.h:988
Role
Definition: groupMerged.h:1132
void addEGPs(const Snippet::Base::EGPVec &egps, const std::string &arrayPrefix, const std::string &varName="")
Definition: groupMerged.h:259
std::vector< DerivedParam > DerivedParamVec
Definition: snippet.h:181
CustomUpdateHostReductionGroupMergedBase(size_t index, const std::string &precision, const BackendBase &backend, const std::vector< std::reference_wrapper< const G >> &groups)
Definition: groupMerged.h:1222
void addScalarField(const std::string &name, GetFieldValueFunc getFieldValue, FieldType fieldType=FieldType::Standard)
Definition: groupMerged.h:220
std::string getKernelSize(const G *group, size_t dimensionIndex, K getKernelSizeFn)
Definition: codeGenUtils.h:167
std::vector< std::string > StringVec
Definition: snippet.h:178
size_t getSize(const std::string &type) const
Get the size of the type.
Definition: backendBase.cc:44
const std::vector< SynapseGroupInternal * > & getSortedArchetypeMergedInSyns() const
Get sorted vectors of merged incoming synapse groups belonging to archetype group.
Definition: groupMerged.h:629
virtual std::string getMergedGroupFieldHostType(const std::string &type) const =0
When generating function calls to push to merged groups, backend without equivalent of Unified Virtua...
VarAccessDuplication
Flags defining how variables should be duplicated across multiple batches.
Definition: varAccess.h:28
void orderNeuronGroupChildren(std::vector< std::vector< T *>> &sortedGroupChildren, G getVectorFunc, H getHashDigestFunc) const
Definition: groupMerged.h:647
void addVars(const Models::Base::VarVec &vars, const std::string &arrayPrefix)
Definition: groupMerged.h:237
size_t getStructArraySize(const BackendBase &backend) const
Definition: groupMerged.h:139
std::string getPostVarIndex(unsigned int batchSize, VarAccessDuplication varDuplication, const std::string &index) const
Definition: groupMerged.h:1092
void addHeterogeneousChildVarInitParams(const Snippet::Base::StringVec &paramNames, const std::vector< std::vector< C >> &sortedGroupChildren, size_t childIndex, size_t varIndex, const std::string &prefix, H isChildParamHeterogeneousFn, V getVarInitialiserFn)
Definition: groupMerged.h:782
std::vector< VarRef > VarRefVec
Definition: models.h:116
void updateHash(H getHashableFn, boost::uuids::detail::sha1 &hash) const
Helper to update hash with the hash of calling getHashableFn on each group.
Definition: groupMerged.h:367
virtual std::string getPointerPrefix() const
Different backends may have different or no pointer prefix (e.g. __global for OpenCL) ...
Definition: backendBase.h:369
std::string getPreWUVarIndex(unsigned int batchSize, VarAccessDuplication varDuplication, const std::string &index) const
Definition: groupMerged.h:1097
bool isParamValueHeterogeneous(size_t index, P getParamValuesFn) const
Helper to test whether parameter values are heterogeneous within merged group.
Definition: groupMerged.h:201
void addVarReferences(const Models::Base::VarRefVec &varReferences, const std::string &arrayPrefix, V getVarRefFn)
Definition: groupMerged.h:246
std::tuple< std::string, std::string, GetFieldValueFunc, FieldType > Field
Definition: groupMerged.h:52
This variable is read-write.
std::string getKernelSize(size_t dimensionIndex) const
Get expression for kernel size in dimension (may be literal or group->kernelSizeXXX) ...
Definition: groupMerged.h:1073
virtual VarRefVec getVarRefs() const
Gets names and types (as strings) of model variable references.
Definition: customUpdateModels.h:38
static const std::string name
Definition: groupMerged.h:1292
void generateRunnerBase(const BackendBase &backend, CodeStream &definitionsInternal, CodeStream &definitionsInternalFunc, CodeStream &definitionsInternalVar, CodeStream &runnerVarDecl, CodeStream &runnerMergedStructAlloc, const std::string &name, bool host=false) const
Definition: groupMerged.h:437
void updateParamHash(R isParamReferencedFn, V getValueFn, boost::uuids::detail::sha1 &hash) const
Definition: groupMerged.h:375
virtual DerivedParamVec getDerivedParams() const
Definition: snippet.h:191
std::string getPreVarIndex(unsigned int batchSize, VarAccessDuplication varDuplication, const std::string &index) const
Definition: groupMerged.h:1087
G GroupInternal
Definition: groupMerged.h:50
virtual std::string getDeviceVarPrefix() const
Definition: backendBase.h:362
bool isKernelSizeHeterogeneous(size_t dimensionIndex) const
Is kernel size heterogeneous in this dimension?
Definition: groupMerged.h:1067
std::string getPostWUVarIndex(unsigned int batchSize, VarAccessDuplication varDuplication, const std::string &index) const
Definition: groupMerged.h:1102
void generateStructFieldArgumentDefinitions(CodeStream &os, const BackendBase &backend) const
Definition: groupMerged.h:126
const std::vector< double > & getDerivedParams() const
Definition: snippet.h:265
std::string getPostISynIndex(unsigned int batchSize, const std::string &index) const
Definition: groupMerged.h:1115
m
Definition: genn_model.py:117
void updateChildVarInitDerivedParamsHash(const std::vector< std::vector< C >> &sortedGroupChildren, size_t childIndex, size_t varIndex, R isChildDerivedParamReferencedFn, V getVarInitialiserFn, boost::uuids::detail::sha1 &hash) const
Definition: groupMerged.h:909
Definition: codeStream.h:94
std::vector< EGP > EGPVec
Definition: snippet.h:179