NEURON
codegen_neuron_cpp_visitor.cpp
Go to the documentation of this file.
1 /*
2  * Copyright 2023 Blue Brain Project, EPFL.
3  * See the top-level LICENSE file for details.
4  *
5  * SPDX-License-Identifier: Apache-2.0
6  */
7 
9 
10 #include <algorithm>
11 #include <chrono>
12 #include <cmath>
13 #include <ctime>
14 #include <optional>
15 #include <regex>
16 #include <stdexcept>
17 
18 #include "ast/all.hpp"
19 #include "ast/procedure_block.hpp"
22 #include "codegen_naming.hpp"
23 #include "config/config.h"
24 #include "parser/c11_driver.hpp"
25 #include "solver/solver.hpp"
26 #include "utils/string_utils.hpp"
30 
31 namespace nmodl {
32 namespace codegen {
33 
34 using namespace ast;
35 
36 using visitor::RenameVisitor;
37 using visitor::VarUsageVisitor;
38 
40 
41 /****************************************************************************************/
42 /* Generic information getters */
43 /****************************************************************************************/
44 
45 
47  return "NEURON";
48 }
49 
50 
52  return "C++ (api-compatibility)";
53 }
54 
55 
57  return "_check_table_thread";
58 }
59 
61  return true;
62 }
63 
64 
65 /****************************************************************************************/
66 /* Common helper routines accross codegen functions */
67 /****************************************************************************************/
68 
70  return get_index_from_name(codegen_float_variables, name);
71 }
72 
73 
74 int CodegenNeuronCppVisitor::position_of_int_var(const std::string& name) const {
75  return get_index_from_name(codegen_int_variables, name);
76 }
77 
78 
79 /****************************************************************************************/
80 /* Backend specific routines */
81 /****************************************************************************************/
82 
84  if (optimize_ionvar_copies) {
85  throw std::runtime_error("Not implemented.");
86  }
87  return false;
88 }
89 
90 
91 /****************************************************************************************/
92 /* Printing routines for code generation */
93 /****************************************************************************************/
94 
95 
97  if (info.point_process) {
98  printer->add_multi_line(R"CODE(
99  /* Point Process specific functions */
100  static void* _hoc_create_pnt(Object* _ho) {
101  return create_point_process(_pointtype, _ho);
102  }
103  )CODE");
104  printer->push_block("static void _hoc_destroy_pnt(void* _vptr)");
105  if (info.is_watch_used() || info.for_netcon_used) {
106  printer->add_line("Prop* _prop = ((Point_process*)_vptr)->prop;");
107  printer->push_block("if (_prop)");
108  printer->add_line("Datum* _ppvar = _nrn_mechanism_access_dparam(_prop);");
109  if (info.is_watch_used()) {
110  printer->fmt_line("_nrn_free_watch(_ppvar, {}, {});",
111  info.watch_count,
112  info.is_watch_used());
113  }
114  if (info.for_netcon_used) {
115  auto fornetcon_data = get_variable_name("fornetcon_data", false);
116  printer->fmt_line("_nrn_free_fornetcon(&{});", fornetcon_data);
117  }
118  printer->pop_block();
119  }
120  printer->add_line("destroy_point_process(_vptr);");
121  printer->pop_block();
122  printer->add_multi_line(R"CODE(
123  static double _hoc_loc_pnt(void* _vptr) {
124  return loc_point_process(_pointtype, _vptr);
125  }
126 
127  static double _hoc_has_loc(void* _vptr) {
128  return has_loc_point(_vptr);
129  }
130 
131  static double _hoc_get_loc_pnt(void* _vptr) {
132  return (get_loc_point_process(_vptr));
133  }
134  )CODE");
135  }
136 }
137 
138 
140  if (info.table_count == 0) {
141  return;
142  }
143 
144  // print declarations of `check_*` functions
145  for (const auto& function: info.functions_with_table) {
146  auto name = function->get_node_name();
147  auto internal_params = internal_method_parameters();
148  printer->fmt_line("void {}({});",
149  table_update_function_name(name),
150  get_parameter_str(internal_params));
151  }
152 
153  ParamVector args = {{"", "Memb_list*", "", "_ml"},
154  {"", "size_t", "", "id"},
155  {"", "Datum*", "", "_ppvar"},
156  {"", "Datum*", "", "_thread"},
157  {"", "double*", "", "_globals"},
158  {"", "NrnThread*", "", "nt"},
159  {"", "int", "", "_type"},
160  {"", "const _nrn_model_sorted_token&", "", "_sorted_token"}};
161 
162  // definition of `_check_table_thread` function
163  // signature must be same as the `nrn_thread_table_check_t` type
164  printer->fmt_line("static void {}({})", table_thread_function_name(), get_parameter_str(args));
165  printer->push_block();
166  printer->add_line("_nrn_mechanism_cache_range _lmc{_sorted_token, *nt, *_ml, _type};");
167  printer->fmt_line("auto inst = make_instance_{}(&_lmc);", info.mod_suffix);
168  if (!info.artificial_cell) {
169  printer->fmt_line("auto node_data = make_node_data_{}(*nt, *_ml);", info.mod_suffix);
170  }
171  if (!codegen_thread_variables.empty()) {
172  printer->fmt_line("auto _thread_vars = {}(_thread[{}].get<double*>());",
173  thread_variables_struct(),
174  info.thread_var_thread_id);
175  }
176 
177  for (const auto& function: info.functions_with_table) {
178  auto method_name = function->get_node_name();
179  auto method_args = get_arg_str(internal_method_parameters());
180  printer->fmt_line("{}({});", table_update_function_name(method_name), method_args);
181  }
182  printer->pop_block();
183 }
184 
185 
187  printer->add_line("/* Neuron setdata functions */");
188  printer->add_line("extern void _nrn_setdata_reg(int, void(*)(Prop*));");
189  printer->push_block("static void _setdata(Prop* _prop)");
190  if (!info.point_process) {
191  printer->add_multi_line(R"CODE(
192  _extcall_prop = _prop;
193  _prop_id = _nrn_get_prop_id(_prop);
194  )CODE");
195  }
196  printer->pop_block();
197 
198  if (info.point_process) {
199  printer->push_block("static void _hoc_setdata(void* _vptr)");
200  printer->add_multi_line(R"CODE(
201  Prop* _prop;
202  _prop = ((Point_process*)_vptr)->prop;
203  _setdata(_prop);
204  )CODE");
205  } else {
206  printer->push_block("static void _hoc_setdata()");
207  printer->add_multi_line(R"CODE(
208  Prop *_prop = hoc_getdata_range(mech_type);
209  _setdata(_prop);
210  hoc_retpushx(1.);
211  )CODE");
212  }
213  printer->pop_block();
214 }
215 
216 
218  printer->add_newline(2);
219 
220  auto print_decl = [this](const auto& callables) {
221  for (const auto& node: callables) {
222  print_function_declaration(*node, node->get_node_name());
223  printer->add_text(';');
224  printer->add_newline();
225  }
226  };
227 
228  printer->add_line("/* Mechanism procedures and functions */");
229  print_decl(info.functions);
230  print_decl(info.procedures);
231 
232  for (const auto& node: info.function_tables) {
233  auto [params, table_params] = function_table_parameters(*node);
234  printer->fmt_line("double {}({});",
235  method_name(node->get_node_name()),
236  get_parameter_str(params));
237  printer->fmt_line("double {}({});",
238  method_name("table_" + node->get_node_name()),
239  get_parameter_str(table_params));
240  }
241 }
242 
243 
245  const ast::Block& node,
246  const std::string& name,
247  const std::unordered_set<CppObjectSpecifier>& specifiers) {
248  printer->add_newline(2);
249  print_function_declaration(node, name, specifiers);
250  printer->add_text(" ");
251  printer->push_block();
252 
253  // function requires return variable declaration
254  if (node.is_function_block()) {
255  auto type = default_float_data_type();
256  printer->fmt_line("{} ret_{} = 0.0;", type, name);
257  } else {
258  printer->fmt_line("int ret_{} = 0;", name);
259  }
260 
261  if (info.mod_suffix != "nothing" && !info.artificial_cell) {
262  printer->add_line(
263  "double v = node_data.node_voltages ? "
264  "node_data.node_voltages[node_data.nodeindices[id]] : 0.0;");
265  }
266 
267  print_statement_block(*node.get_statement_block(), false, false);
268  printer->fmt_line("return ret_{};", name);
269  printer->pop_block();
270 }
271 
272 
274  auto name = node.get_node_name();
275  if (info.function_uses_table(name)) {
276  auto new_name = "f_" + name;
277  print_function_or_procedure(node,
278  new_name,
280  print_table_check_function(node);
281  print_table_replacement_function(node);
282  } else {
283  print_function_or_procedure(node, name);
284  }
285 }
286 
288  const ast::Block* function_or_procedure_block,
289  InterpreterWrapper wrapper_type) {
290  const auto block_name = function_or_procedure_block->get_node_name();
291 
292  const auto get_func_call_str = [&]() {
293  const auto& params = function_or_procedure_block->get_parameters();
294  const auto func_proc_name = block_name + "_" + info.mod_suffix;
295  std::vector<std::string> args;
296  args.reserve(params.size());
297  for (int i = 0; i < params.size(); ++i) {
298  args.push_back(fmt::format("*getarg({})", i + 1));
299  }
300 
301  auto internal_args = internal_method_arguments();
302  return fmt::format("{}({})",
303  func_proc_name,
304  stringutils::join_arguments(internal_args,
305  fmt::format("{}", fmt::join(args, ", "))));
306  };
307 
308  printer->add_line("double _r = 0.0;");
309  if (function_or_procedure_block->is_function_block()) {
310  printer->add_indent();
311  printer->fmt_text("_r = {};", get_func_call_str());
312  printer->add_newline();
313  } else {
314  printer->add_line("_r = 1.;");
315  printer->fmt_line("{};", get_func_call_str());
316  }
317  if (info.point_process || wrapper_type != InterpreterWrapper::HOC) {
318  printer->add_line("return(_r);");
319  } else if (wrapper_type == InterpreterWrapper::HOC) {
320  printer->add_line("hoc_retpushx(_r);");
321  }
322 }
323 
325  const ast::Block* function_or_procedure_block,
326  InterpreterWrapper wrapper_type) {
327  if (info.mod_suffix == "nothing") {
328  return;
329  }
330 
331  const auto block_name = function_or_procedure_block->get_node_name();
332  printer->add_multi_line(R"CODE(
333  Datum* _ppvar;
334  Datum* _thread;
335  NrnThread* nt;
336  )CODE");
337 
338  std::string prop_name;
339  if (info.point_process) {
340  printer->add_multi_line(R"CODE(
341  auto* const _pnt = static_cast<Point_process*>(_vptr);
342  auto* const _p = _pnt->prop;
343  if (!_p) {
344  hoc_execerror("POINT_PROCESS data instance not valid", nullptr);
345  }
346  _nrn_mechanism_cache_instance _lmc{_p};
347  size_t const id{};
348  _ppvar = _nrn_mechanism_access_dparam(_p);
349  _thread = _extcall_thread.data();
350  nt = static_cast<NrnThread*>(_pnt->_vnt);
351  )CODE");
352 
353  prop_name = "_p";
354  } else if (wrapper_type == InterpreterWrapper::HOC) {
355  if (program_symtab->lookup(block_name)->has_all_properties(NmodlType::use_range_ptr_var)) {
356  printer->push_block("if (!_prop_id)");
357  printer->fmt_line(
358  "hoc_execerror(\"No data for {}_{}. Requires prior call to setdata_{} and that the "
359  "specified mechanism instance still be in existence.\", nullptr);",
360  function_or_procedure_block->get_node_name(),
361  info.mod_suffix,
362  info.mod_suffix);
363  printer->pop_block();
364  printer->add_line("Prop* _local_prop = _extcall_prop;");
365  } else {
366  printer->add_line("Prop* _local_prop = _prop_id ? _extcall_prop : nullptr;");
367  }
368  printer->add_multi_line(R"CODE(
369  _nrn_mechanism_cache_instance _lmc{_local_prop};
370  size_t const id{};
371  _ppvar = _local_prop ? _nrn_mechanism_access_dparam(_local_prop) : nullptr;
372  _thread = _extcall_thread.data();
373  nt = nrn_threads;
374  )CODE");
375  prop_name = "_local_prop";
376  } else { // wrapper_type == InterpreterWrapper::Python
377  printer->add_multi_line(R"CODE(
378  _nrn_mechanism_cache_instance _lmc{_prop};
379  size_t const id = 0;
380  _ppvar = _nrn_mechanism_access_dparam(_prop);
381  _thread = _extcall_thread.data();
382  nt = nrn_threads;
383  )CODE");
384  prop_name = "_prop";
385  }
386 
387  printer->fmt_line("auto inst = make_instance_{}({} ? &_lmc : nullptr);",
388  info.mod_suffix,
389  prop_name);
390  if (!info.artificial_cell) {
391  printer->fmt_line("auto node_data = make_node_data_{}({});", info.mod_suffix, prop_name);
392  }
393  if (!codegen_thread_variables.empty()) {
394  printer->fmt_line("auto _thread_vars = {}(_thread[{}].get<double*>());",
395  thread_variables_struct(),
396  info.thread_var_thread_id);
397  }
398  if (info.function_uses_table(block_name)) {
399  printer->fmt_line("{}({});",
400  table_update_function_name(block_name),
401  internal_method_arguments());
402  }
403 }
404 
405 
407  const ast::Block* function_or_procedure_block,
408  InterpreterWrapper wrapper_type) {
409  const auto block_name = function_or_procedure_block->get_node_name();
410  if (wrapper_type == InterpreterWrapper::HOC) {
411  return hoc_function_signature(block_name);
412  } else {
413  return py_function_signature(block_name);
414  }
415 }
416 
417 void CodegenNeuronCppVisitor::print_hoc_py_wrapper(const ast::Block* function_or_procedure_block,
418  InterpreterWrapper wrapper_type) {
419  if (info.point_process && wrapper_type == InterpreterWrapper::Python) {
420  return;
421  }
422 
423  printer->push_block(hoc_py_wrapper_signature(function_or_procedure_block, wrapper_type));
424 
425  print_hoc_py_wrapper_setup(function_or_procedure_block, wrapper_type);
426  print_hoc_py_wrapper_call_impl(function_or_procedure_block, wrapper_type);
427 
428  printer->pop_block();
429 }
430 
431 
433  auto print_wrappers = [this](const auto& callables) {
434  for (const auto& callable: callables) {
435  print_hoc_py_wrapper(callable, InterpreterWrapper::HOC);
436  print_hoc_py_wrapper(callable, InterpreterWrapper::Python);
437  }
438  };
439 
440  print_wrappers(info.procedures);
441  print_wrappers(info.functions);
442 
443  for (const auto& node: info.function_tables) {
444  auto name = node->get_node_name();
445  auto table_name = "table_" + node->get_node_name();
446 
447 
448  auto args = std::vector<std::string>{};
449  for (size_t i = 0; i < node->get_parameters().size(); ++i) {
450  args.push_back(fmt::format("*getarg({})", i + 1));
451  }
452 
453  // HOC
454  std::string return_statement = info.point_process ? "return _ret;" : "hoc_retpushx(_ret);";
455 
456  printer->fmt_push_block("{}", hoc_function_signature(name));
457  printer->fmt_line("double _ret = {}({});", method_name(name), fmt::join(args, ", "));
458  printer->add_line(return_statement);
459  printer->pop_block();
460 
461  printer->fmt_push_block("{}", hoc_function_signature(table_name));
462  printer->fmt_line("double _ret = {}();", method_name(table_name));
463  printer->add_line(return_statement);
464  printer->pop_block();
465 
466  // Python
467  printer->fmt_push_block("{}", py_function_signature(name));
468  printer->fmt_line("return {}({});", method_name(name), fmt::join(args, ", "));
469  printer->pop_block();
470 
471  printer->fmt_push_block("{}", py_function_signature(table_name));
472  printer->fmt_line("return {}();", method_name(table_name));
473  printer->pop_block();
474  }
475 }
476 
478  return ParamVector{{"", "ldifusfunc2_t", "", "_f"},
479  {"", "const _nrn_model_sorted_token&", "", "_sorted_token"},
480  {"", "NrnThread&", "", "_nt"}};
481 }
482 
483 
485  return ParamVector{{"", "int", "", "_i"},
486  {"", "Memb_list*", "", "_ml_arg"},
487  {"", "size_t", "", "id"},
488  {"", "Datum*", "", "_ppvar"},
489  {"", "double*", "", "_pdvol"},
490  {"", "double*", "", "_pdfcdc"},
491  {"", "Datum*", "", "/* _thread */"},
492  {"", "NrnThread*", "", "nt"},
493  {"", "const _nrn_model_sorted_token&", "", "_sorted_token"}};
494 }
495 
497  auto coeff_callback_name = [](const std::string& var_name) {
498  return fmt::format("_diffusion_coefficient_{}", var_name);
499  };
500 
501  auto space_name = [](const std::string& var_name) {
502  return fmt::format("_diffusion_space_{}", var_name);
503  };
504 
505  for (auto [var_name, values]: info.longitudinal_diffusion_info) {
506  printer->fmt_line("static void* {};", space_name(var_name));
507  printer->fmt_push_block("static double {}({})",
508  coeff_callback_name(var_name),
509  get_parameter_str(ldifusfunc3_parameters()));
510 
511  print_entrypoint_setup_code_from_memb_list();
512 
513  auto volume_expr = values.volume("_i");
514  auto mu_expr = values.diffusion_rate("_i");
515 
516  printer->add_indent();
517  printer->add_text("*_pdvol= ");
518  volume_expr->accept(*this);
519  printer->add_text(";");
520  printer->add_newline();
521 
522  printer->add_line("*_pdfcdc = 0.0;");
523  printer->add_indent();
524  printer->add_text("return ");
525  mu_expr->accept(*this);
526  printer->add_text(";");
527  printer->add_newline();
528 
529  printer->pop_block();
530  }
531 
532  printer->fmt_push_block("static void _apply_diffusion_function({})",
533  get_parameter_str(ldifusfunc1_parameters()));
534  for (auto [var_name, values]: info.longitudinal_diffusion_info) {
535  auto var = program_symtab->lookup(var_name);
536  size_t array_size = var->get_length();
537  printer->fmt_push_block("for(size_t _i = 0; _i < {}; ++_i)", array_size);
538  printer->fmt_line(
539  "(*_f)(mech_type, {}, &{}, _i, /* x pos */ {}, /* Dx pos */ {}, _sorted_token, _nt);",
540  coeff_callback_name(var_name),
541  space_name(var_name),
542  position_of_float_var(var_name),
543  position_of_float_var("D" + var_name));
544  printer->pop_block();
545  }
546  printer->pop_block();
547  printer->add_newline();
548 }
549 
550 /****************************************************************************************/
551 /* Code-specific helper routines */
552 /****************************************************************************************/
553 
554 void CodegenNeuronCppVisitor::add_variable_tqitem(std::vector<IndexVariableInfo>& variables) {
555  if (info.net_send_used) {
556  variables.emplace_back(make_symbol(naming::TQITEM_VARIABLE), false, false, true);
557  variables.back().is_constant = true;
558  info.tqitem_index = static_cast<int>(variables.size() - 1);
559  }
560 }
561 
563  std::vector<IndexVariableInfo>& variables) {
564  variables.emplace_back(make_symbol(naming::POINT_PROCESS_VARIABLE), false, false, true);
565  variables.back().is_constant = true;
566 }
567 
569  const auto& args = internal_method_parameters();
570  return get_arg_str(args);
571 }
572 
573 
575  if (info.mod_suffix == "nothing") {
576  return {};
577  }
578 
579  ParamVector params;
580  params.emplace_back("", "_nrn_mechanism_cache_range&", "", "_lmc");
581  params.emplace_back("", fmt::format("{}&", instance_struct()), "", "inst");
582  if (!info.artificial_cell) {
583  params.emplace_back("", fmt::format("{}&", node_data_struct()), "", "node_data");
584  }
585  params.emplace_back("", "size_t", "", "id");
586  params.emplace_back("", "Datum*", "", "_ppvar");
587  params.emplace_back("", "Datum*", "", "_thread");
588  if (!codegen_thread_variables.empty()) {
589  params.emplace_back("", fmt::format("{}&", thread_variables_struct()), "", "_thread_vars");
590  }
591  params.emplace_back("", "NrnThread*", "", "nt");
592  return params;
593 }
594 
595 
596 /// TODO: Edit for NEURON
598  return {};
599 }
600 
601 
602 /// TODO: Edit for NEURON
604  bool table) noexcept {
605  return {};
606 }
607 
608 
610  return internal_method_parameters();
611 }
612 
613 
615  return {{"", "Memb_list*", "", "_ml"},
616  {"", "size_t", "", "_iml"},
617  {"", "Datum*", "", "_ppvar"},
618  {"", "Datum*", "", "_thread"},
619  {"", "double*", "", "_globals"},
620  {"", "NrnThread*", "", "_nt"}};
621 }
622 
623 
624 /// TODO: Edit for NEURON
626  return {};
627 }
628 
629 
630 /// TODO: Edit for NEURON
632  return {};
633 }
634 
635 std::pair<CodegenNeuronCppVisitor::ParamVector, CodegenNeuronCppVisitor::ParamVector>
637  auto params = ParamVector{};
638 
639  for (const auto& i: node.get_parameters()) {
640  params.emplace_back("", "double", "", i->get_node_name());
641  }
642  return {params, {}};
643 }
644 
645 
646 /** Map of the non-(global/top-local) LOCAL variables.
647  *
648  * The map associates the name in the MOD file, e.g. `a` with
649  * the current name of that LOCAL variable, e.g. `a_r_4`.
650  *
651  * auto map = get_nonglobal_local_variable_names();
652  * assert map["a"] == "a_r_4";
653  *
654  * The two names can differ, because an early pass makes all
655  * names unique by renaming local variables.
656  */
657 std::unordered_map<std::string, std::string> get_nonglobal_local_variable_names(
658  const symtab::SymbolTable& symtab) {
659  if (symtab.global_scope()) {
660  return {};
661  }
662 
663  auto local_variables = symtab.get_variables(NmodlType::local_var);
664  auto parent_symtab = symtab.get_parent_table();
665  if (parent_symtab == nullptr) {
666  throw std::runtime_error(
667  "Internal NMODL error: non top-level symbol table doesn't have a parent.");
668  }
669 
670  auto variable_names = get_nonglobal_local_variable_names(*parent_symtab);
671 
672  for (const auto& symbol: local_variables) {
673  auto status = symbol->get_status();
674  bool is_renamed = (status & symtab::syminfo::Status::renamed) !=
676  auto current_name = symbol->get_name();
677  auto mod_name = is_renamed ? symbol->get_original_name() : current_name;
678 
679  variable_names[mod_name] = current_name;
680  }
681 
682  return variable_names;
683 }
684 
685 
687  const ast::Verbatim& node,
688  const std::string& verbatim) {
689  // Note, the logic for reducing the number of macros printed, aims to
690  // improve legibility of the generated code by reducing number of lines of
691  // code. It would be correct to print all macros, because that's
692  // essentially what NOCMODL does. Therefore, the logic isn't sharp (and
693  // doesn't have to be).
694 
695  std::vector<std::string> macros_defined;
696  auto print_macro = [this, &verbatim, &macros_defined](const std::string& macro_name,
697  const std::string& macro_value) {
698  if (verbatim.find(macro_name) != std::string::npos) {
699  printer->fmt_line("#define {} {}", macro_name, macro_value);
700  macros_defined.push_back(macro_name);
701  }
702  };
703 
704  printer->add_line("// Setup for VERBATIM");
705  for (const auto& var: codegen_float_variables) {
706  auto name = get_name(var);
707  print_macro(name, get_variable_name(name));
708  }
709 
710  for (const auto& var: codegen_int_variables) {
711  auto name = get_name(var);
712  std::string macro_value = get_variable_name(name);
713  print_macro(name, macro_value);
714  if (verbatim.find("_p_" + name) != std::string::npos) {
715  print_macro("_p_" + name, get_pointer_name(name));
716  }
717  }
718 
719  for (const auto& var: codegen_global_variables) {
720  auto name = get_name(var);
721  print_macro(name, get_variable_name(name));
722  }
723 
724  for (const auto& var: codegen_thread_variables) {
725  auto name = get_name(var);
726  print_macro(name, get_variable_name(name));
727  }
728 
729  for (const auto& func: info.functions) {
730  auto name = get_name(func);
731  print_macro(name, method_name(name));
732  print_macro(fmt::format("_l{}", name), fmt::format("ret_{}", name));
733  }
734 
735  for (const auto& proc: info.procedures) {
736  auto name = get_name(proc);
737  print_macro(name, method_name(name));
738  }
739 
740 
741  auto symtab = node.get_parent()->get_symbol_table();
742  auto locals = get_nonglobal_local_variable_names(*symtab);
743  for (const auto& [mod_name, current_name]: locals) {
744  print_macro(fmt::format("_l{}", mod_name), get_variable_name(current_name));
745  }
746 
747  print_macro("v", get_variable_name("v"));
748  print_macro(naming::NTHREAD_T_VARIABLE, "nt->_t");
749  print_macro(naming::NTHREAD_DT_VARIABLE, "nt->_dt");
750  print_macro("_nt", "nt");
751  print_macro("_tqitem", "tqitem");
752 
753  auto print_args_macro = [this, print_macro](const std::string& macro_basename,
754  const ParamVector& params) {
755  print_macro("_" + macro_basename + "_", get_arg_str(params));
756  print_macro("_" + macro_basename + "comma_", get_arg_str(params) + ",");
757  print_macro("_" + macro_basename + "proto_", get_parameter_str(params));
758  print_macro("_" + macro_basename + "protocomma_", get_parameter_str(params) + ",");
759  };
760 
761  print_args_macro("internalthreadargs", internalthreadargs_parameters());
762  print_args_macro("threadargs", threadargs_parameters());
763 
764  return macros_defined;
765 }
766 
767 
769  const std::vector<std::string>& macros_defined) {
770  for (const auto& macro: macros_defined) {
771  printer->fmt_line("#undef {}", macro);
772  }
773  printer->add_line("// End of cleanup for VERBATIM");
774 }
775 
776 
777 std::string CodegenNeuronCppVisitor::process_verbatim_text(const std::string& verbatim) {
779  driver.scan_string(verbatim);
780  auto tokens = driver.all_tokens();
781  std::string result;
782  for (size_t i = 0; i < tokens.size(); i++) {
783  auto token = tokens[i];
784 
785  // check if we have function call in the verbatim block where
786  // function is defined in the same mod file
787  if (program_symtab->is_method_defined(token) && tokens[i + 1] == "(") {
788  result += token + "(";
789  if (tokens[i + 2] == naming::THREAD_ARGS) {
791  } else if (tokens[i + 2] == naming::THREAD_ARGS_COMMA) {
793  } else {
794  result += tokens[i + 2];
795  }
796 
797  i += 2;
798  } else {
799  result += token;
800  }
801  }
802  return result;
803 }
804 
805 
807  const auto& verbatim_code = node.get_statement()->eval();
808  auto massaged_verbatim = process_verbatim_text(verbatim_code);
809 
810  auto macros_defined = print_verbatim_setup(node, massaged_verbatim);
811  printer->add_line("// Begin VERBATIM");
812  const auto& lines = stringutils::split_string(massaged_verbatim, '\n');
813  for (const auto& line: lines) {
814  printer->add_line(line);
815  }
816  printer->add_line("// End VERBATIM");
817  print_verbatim_cleanup(macros_defined);
818 }
819 
820 
821 /// TODO: Write for NEURON
823  return {};
824 };
825 
826 
828  const std::string& function_or_procedure_name) const {
829  return fmt::format("_hoc_{}", function_or_procedure_name);
830 }
831 
832 
834  const std::string& function_or_procedure_name) const {
835  return fmt::format("static {} {}({})",
836  info.point_process ? "double" : "void",
837  hoc_function_name(function_or_procedure_name),
838  info.point_process ? "void * _vptr" : "");
839 }
840 
841 
843  const std::string& function_or_procedure_name) const {
844  return fmt::format("_npy_{}", function_or_procedure_name);
845 }
846 
847 
849  const std::string& function_or_procedure_name) const {
850  return fmt::format("static double {}(Prop* _prop)",
851  py_function_name(function_or_procedure_name));
852 }
853 
854 
855 /****************************************************************************************/
856 /* Code-specific printing routines for code generation */
857 /****************************************************************************************/
858 
860  return "neuron";
861 }
862 
863 
865  std::vector<ShadowUseStatement>& statements,
866  const Ion& ion,
867  const std::string& /* concentration */) {
868  auto ion_name = ion.name;
869  int dparam_index = get_int_variable_index(fmt::format("style_{}", ion_name));
870 
871  auto style_name = fmt::format("_style_{}", ion_name);
872  auto style_stmt = fmt::format("int {} = *(_ppvar[{}].get<int*>())", style_name, dparam_index);
873  statements.push_back(ShadowUseStatement{style_stmt, "", ""});
874 
875 
876  auto wrote_conc_stmt = fmt::format("nrn_wrote_conc(_{}_sym, {}, {}, {}, {})",
877  ion_name,
878  get_variable_name(ion.rev_potential_pointer_name()),
879  get_variable_name(ion.intra_conc_pointer_name()),
880  get_variable_name(ion.extra_conc_pointer_name()),
881  style_name);
882  statements.push_back(ShadowUseStatement{wrote_conc_stmt, "", ""});
883 }
884 
885 /****************************************************************************************/
886 /* Routines for returning variable name */
887 /****************************************************************************************/
888 
889 
891  bool use_instance) const {
892  if (!use_instance) {
893  throw std::runtime_error("Printing non-instance variables is not implemented.");
894  }
895 
896  auto name = symbol->get_name();
897  auto dimension = symbol->get_length();
898  if (symbol->is_array()) {
899  return fmt::format("(inst.{}+id*{})", name, dimension);
900  } else {
901  return fmt::format("inst.{}[id]", name);
902  }
903 }
904 
905 
907  const std::string& name,
908  bool use_instance) const {
909  auto position = position_of_int_var(name);
910 
911  if (info.semantics[position].name == naming::RANDOM_SEMANTIC) {
912  return fmt::format("(nrnran123_State*) _ppvar[{}].literal_value<void*>()", position);
913  }
914 
915  if (info.semantics[position].name == naming::FOR_NETCON_SEMANTIC) {
916  return fmt::format("_ppvar[{}].literal_value<void*>()", position);
917  }
918 
919  if (info.semantics[position].name == naming::POINTER_SEMANTIC) {
920  return fmt::format("(*_ppvar[{}].get<double*>())", position);
921  }
922 
923  if (symbol.is_index) {
924  if (use_instance) {
925  throw std::runtime_error("Not implemented. [wiejo]");
926  // return fmt::format("inst->{}[{}]", name, position);
927  }
928  throw std::runtime_error("Not implemented. [ncuwi]");
929  // return fmt::format("indexes[{}]", position);
930  }
931  if (symbol.is_integer) {
932  if (use_instance) {
933  return fmt::format("inst.{}[id]", name);
934  }
935  return fmt::format("_ppvar[{}]", position);
936  }
937  if (use_instance) {
938  return fmt::format("(*inst.{}[id])", name);
939  }
940 
941 
942  throw std::runtime_error("Not implemented. [nvueir]");
943 }
944 
945 
947  bool use_instance) const {
948  auto i_var = var_info.offset;
949  auto var_name = var_info.symbol->get_name();
950 
951  if (use_instance) {
952  if (var_info.symbol->is_array()) {
953  return fmt::format("(_thread_vars.{}_ptr(id))", var_name);
954  } else {
955  return fmt::format("_thread_vars.{}(id)", var_name);
956  }
957  } else {
958  if (var_info.symbol->is_array()) {
959  return fmt::format("({}.thread_data + {})", global_struct_instance(), i_var);
960  } else {
961  return fmt::format("{}.thread_data[{}]", global_struct_instance(), i_var);
962  }
963  }
964 }
965 
966 
968  bool use_instance) const {
969  if (use_instance) {
970  return fmt::format("inst.{}->{}", naming::INST_GLOBAL_MEMBER, symbol->get_name());
971  } else {
972  return fmt::format("{}.{}", global_struct_instance(), symbol->get_name());
973  }
974 }
975 
976 
977 std::string CodegenNeuronCppVisitor::get_pointer_name(const std::string& name) const {
978  auto name_comparator = [&name](const auto& sym) { return name == get_name(sym); };
979 
980  auto var =
981  std::find_if(codegen_int_variables.begin(), codegen_int_variables.end(), name_comparator);
982 
983  if (var == codegen_int_variables.end()) {
984  throw std::runtime_error("Only integer variables have a 'pointer name'.");
985  }
986  auto position = position_of_int_var(name);
987  return fmt::format("_ppvar[{}].literal_value<void*>()", position);
988 }
989 
990 
991 std::string CodegenNeuronCppVisitor::get_variable_name(const std::string& name,
992  bool use_instance) const {
993  std::string varname = update_if_ion_variable_name(name);
994  if (!info.artificial_cell && varname == "v") {
996  }
997 
998  auto name_comparator = [&varname](const auto& sym) { return varname == get_name(sym); };
999 
1001  if (printing_net_receive) {
1002  // In net_receive blocks, the point process is passed in as an
1003  // argument called:
1004  return "_pnt";
1005  }
1006  // The "integer variable" branch will pick up the correct `_ppvar` when
1007  // not printing a NET_RECEIVE block.
1008  }
1009 
1010  if (stringutils::starts_with(name, "_ptable_")) {
1011  return fmt::format("{}.{}", global_struct_instance(), name);
1012  }
1013 
1014  // float variable
1015  auto f = std::find_if(codegen_float_variables.begin(),
1016  codegen_float_variables.end(),
1017  name_comparator);
1018  if (f != codegen_float_variables.end()) {
1019  return float_variable_name(*f, use_instance);
1020  }
1021 
1022  // integer variable
1023  auto i =
1024  std::find_if(codegen_int_variables.begin(), codegen_int_variables.end(), name_comparator);
1025  if (i != codegen_int_variables.end()) {
1026  return int_variable_name(*i, varname, use_instance);
1027  }
1028 
1029  // thread variable
1030  auto t = std::find_if(codegen_thread_variables.begin(),
1031  codegen_thread_variables.end(),
1032  name_comparator);
1033  if (t != codegen_thread_variables.end()) {
1034  return thread_variable_name(*t, use_instance);
1035  }
1036 
1037  // global variable
1038  auto g = std::find_if(codegen_global_variables.begin(),
1039  codegen_global_variables.end(),
1040  name_comparator);
1041  if (g != codegen_global_variables.end()) {
1042  return global_variable_name(*g, use_instance);
1043  }
1044 
1045  if (varname == naming::NTHREAD_DT_VARIABLE) {
1046  return std::string("nt->_") + naming::NTHREAD_DT_VARIABLE;
1047  }
1048 
1049  if (varname == naming::NTHREAD_T_VARIABLE) {
1050  return std::string("nt->_") + naming::NTHREAD_T_VARIABLE;
1051  }
1052 
1053  // external variable
1054  auto e = std::find_if(info.external_variables.begin(),
1055  info.external_variables.end(),
1056  name_comparator);
1057  if (e != info.external_variables.end()) {
1058  return fmt::format("{}()", varname);
1059  }
1060 
1061  auto const iter =
1062  std::find_if(info.neuron_global_variables.begin(),
1063  info.neuron_global_variables.end(),
1064  [&varname](auto const& entry) { return entry.first->get_name() == varname; });
1065  if (iter != info.neuron_global_variables.end()) {
1066  std::string ret;
1067  if (use_instance) {
1068  ret = "*(inst.";
1069  }
1070  ret.append(varname);
1071  if (use_instance) {
1072  ret.append(")");
1073  }
1074  return ret;
1075  }
1076 
1077  // otherwise return original name
1078  return varname;
1079 }
1080 
1081 
1082 /****************************************************************************************/
1083 /* Main printing routines for code generation */
1084 /****************************************************************************************/
1085 
1086 
1088  printer->add_newline();
1089  printer->add_multi_line(R"CODE(
1090  #include <Eigen/Dense>
1091  #include <Eigen/LU>
1092  #include <math.h>
1093  #include <stdio.h>
1094  #include <stdlib.h>
1095  #include <vector>
1096  )CODE");
1097  if (info.eigen_newton_solver_exist) {
1098  printer->add_multi_line(nmodl::solvers::newton_hpp);
1099  }
1100  if (!info.matexp_blocks.empty()) {
1101  printer->add_line("#include <unsupported/Eigen/MatrixFunctions>");
1102  }
1103 }
1104 
1105 
1107  printer->add_newline();
1108  printer->add_multi_line(R"CODE(
1109  #include "mech_api.h"
1110  #include "neuron/cache/mechanism_range.hpp"
1111  #include "nmodlmutex.h"
1112  #include "nrniv_mf.h"
1113  #include "section_fwd.hpp"
1114  )CODE");
1115 }
1116 
1117 
1118 void CodegenNeuronCppVisitor::print_sdlists_init(bool /* print_initializers */) {
1119  /// _initlists() should only be called once by the mechanism registration function
1120  /// (_<mod_file>_reg())
1121  printer->add_newline(2);
1122  printer->push_block("static void _initlists()");
1123  for (auto i = 0; i < info.prime_variables_by_order.size(); ++i) {
1124  const auto& prime_var = info.prime_variables_by_order[i];
1125  /// TODO: Something similar needs to happen for slist/dlist2 but I don't know their usage at
1126  // the moment
1127  /// TODO: We have to do checks and add errors similar to nocmodl in the
1128  // SemanticAnalysisVisitor
1129  if (prime_var->is_array()) {
1130  /// TODO: Needs a for loop here. Look at
1131  // https://github.com/neuronsimulator/nrn/blob/df001a436bcb4e23d698afe66c2a513819a6bfe8/src/nmodl/deriv.cpp#L524
1132  /// TODO: Also needs a test
1133  printer->fmt_push_block("for (int _i = 0; _i < {}; ++_i)", prime_var->get_length());
1134  printer->fmt_line("/* {}[{}] */", prime_var->get_name(), prime_var->get_length());
1135  printer->fmt_line("_slist1[{}+_i] = {{{}, _i}};",
1136  i,
1137  position_of_float_var(prime_var->get_name()));
1138  const auto prime_var_deriv_name = "D" + prime_var->get_name();
1139  printer->fmt_line("/* {}[{}] */", prime_var_deriv_name, prime_var->get_length());
1140  printer->fmt_line("_dlist1[{}+_i] = {{{}, _i}};",
1141  i,
1142  position_of_float_var(prime_var_deriv_name));
1143  printer->pop_block();
1144  } else {
1145  printer->fmt_line("/* {} */", prime_var->get_name());
1146  printer->fmt_line("_slist1[{}] = {{{}, 0}};",
1147  i,
1148  position_of_float_var(prime_var->get_name()));
1149  const auto prime_var_deriv_name = "D" + prime_var->get_name();
1150  printer->fmt_line("/* {} */", prime_var_deriv_name);
1151  printer->fmt_line("_dlist1[{}] = {{{}, 0}};",
1152  i,
1153  position_of_float_var(prime_var_deriv_name));
1154  }
1155  }
1156  printer->pop_block();
1157 }
1158 
1160  auto params = internal_method_parameters();
1161 
1162  return params;
1163 }
1164 
1166  const auto value_initialize = print_initializers ? "{}" : "";
1167 
1168  /// TODO: Print only global variables printed in NEURON
1169  printer->add_newline(2);
1170  printer->add_line("/* NEURON global variables */");
1171  if (info.primes_size != 0) {
1172  printer->fmt_line("static neuron::container::field_index _slist1[{0}], _dlist1[{0}];",
1173  info.primes_size);
1174  }
1175 
1176  for (const auto& ion: info.ions) {
1177  printer->fmt_line("static Symbol* _{}_sym;", ion.name);
1178  }
1179 
1180  if (info.emit_cvode) {
1181  printer->add_line("static Symbol** _atollist;");
1182  printer->push_block("static HocStateTolerance _hoc_state_tol[] =");
1183  printer->add_line("{0, 0}");
1184  printer->pop_block(";");
1185  }
1186 
1187  printer->add_line("static int mech_type;");
1188 
1189  if (info.point_process) {
1190  printer->add_line("static int _pointtype;");
1191  } else {
1192  printer->add_multi_line(R"CODE(
1193  static Prop* _extcall_prop;
1194  /* _prop_id kind of shadows _extcall_prop to allow validity checking. */
1195  static _nrn_non_owning_id_without_container _prop_id{};)CODE");
1196  }
1197 
1198  printer->add_line("static _nrn_mechanism_std_vector<Datum> _extcall_thread;");
1199 
1200  // Start printing the CNRN-style global variables.
1201  auto float_type = default_float_data_type();
1202  printer->add_newline(2);
1203  printer->add_line("/** all global variables */");
1204  printer->fmt_push_block("struct {}", global_struct());
1205 
1206  if (!info.ions.empty()) {
1207  // TODO implement these when needed.
1208  }
1209 
1210  if (!info.thread_variables.empty()) {
1211  size_t prefix_sum = 0;
1212  for (size_t i = 0; i < info.thread_variables.size(); ++i) {
1213  const auto& var = info.thread_variables[i];
1214  codegen_thread_variables.push_back({var, i, prefix_sum});
1215 
1216  prefix_sum += var->get_length();
1217  }
1218  }
1219 
1220 
1221  for (const auto& var: info.global_variables) {
1222  codegen_global_variables.push_back(var);
1223  }
1224 
1225  if (info.vectorize && !info.top_local_variables.empty()) {
1226  size_t prefix_sum = info.thread_var_data_size;
1227  size_t n_thread_vars = codegen_thread_variables.size();
1228  for (size_t i = 0; i < info.top_local_variables.size(); ++i) {
1229  const auto& var = info.top_local_variables[i];
1230  codegen_thread_variables.push_back({var, n_thread_vars + i, prefix_sum});
1231 
1232  prefix_sum += var->get_length();
1233  }
1234  }
1235 
1236  if (!info.vectorize && !info.top_local_variables.empty()) {
1237  for (size_t i = 0; i < info.top_local_variables.size(); ++i) {
1238  const auto& var = info.top_local_variables[i];
1239  codegen_global_variables.push_back(var);
1240  }
1241  }
1242 
1243 
1244  if (!codegen_thread_variables.empty()) {
1245  if (!info.vectorize) {
1246  // MOD files that aren't "VECTORIZED" don't have thread data.
1247  throw std::runtime_error("Found thread variables with `vectorize == false`.");
1248  }
1249 
1250  codegen_global_variables.push_back(make_symbol("thread_data_in_use"));
1251 
1252  auto symbol = make_symbol("thread_data");
1253  auto thread_data_size = info.thread_var_data_size + info.top_local_thread_size;
1254  symbol->set_as_array(thread_data_size);
1255  codegen_global_variables.push_back(symbol);
1256  }
1257 
1258  for (const auto& var: info.state_vars) {
1259  auto name = var->get_name() + "0";
1260  auto symbol = program_symtab->lookup(name);
1261  if (symbol == nullptr) {
1262  codegen_global_variables.push_back(make_symbol(name));
1263  }
1264  }
1265 
1266  for (const auto& var: info.constant_variables) {
1267  codegen_global_variables.push_back(var);
1268  }
1269 
1270  for (const auto& var: codegen_global_variables) {
1271  auto name = var->get_name();
1272  auto length = var->get_length();
1273  if (var->is_array()) {
1274  printer->fmt_line("{} {}[{}] /* TODO init const-array */;", float_type, name, length);
1275  } else {
1276  double value{};
1277  if (auto const& value_ptr = var->get_value()) {
1278  value = *value_ptr;
1279  }
1280  printer->fmt_line("{} {}{};",
1281  float_type,
1282  name,
1283  print_initializers ? fmt::format("{{{:g}}}", value) : std::string{});
1284  }
1285  }
1286 
1287  if (info.table_count > 0) {
1288  // basically the same code as coreNEURON uses
1289  printer->fmt_line("double usetable{};", print_initializers ? "{1}" : "");
1290  codegen_global_variables.push_back(make_symbol(naming::USE_TABLE_VARIABLE));
1291 
1292  for (const auto& block: info.functions_with_table) {
1293  const auto& name = block->get_node_name();
1294  printer->fmt_line("{} tmin_{}{};", float_type, name, value_initialize);
1295  printer->fmt_line("{} mfac_{}{};", float_type, name, value_initialize);
1296  codegen_global_variables.push_back(make_symbol("tmin_" + name));
1297  codegen_global_variables.push_back(make_symbol("mfac_" + name));
1298  }
1299 
1300  for (const auto& variable: info.table_statement_variables) {
1301  auto const name = "t_" + variable->get_name();
1302  auto const num_values = variable->get_num_values();
1303  if (variable->is_array()) {
1304  int array_len = variable->get_length();
1305  printer->fmt_line(
1306  "{} {}[{}][{}]{};", float_type, name, array_len, num_values, value_initialize);
1307  } else {
1308  printer->fmt_line("{} {}[{}]{};", float_type, name, num_values, value_initialize);
1309  }
1310  codegen_global_variables.push_back(make_symbol(name));
1311  }
1312  }
1313 
1314  print_global_struct_function_table_ptrs();
1315 
1316  if (info.vectorize && info.thread_data_index) {
1317  // TODO compare CoreNEURON something extcall stuff.
1318  // throw std::runtime_error("Not implemented, global vectorize something else.");
1319  }
1320 
1321  if (info.diam_used) {
1322  printer->fmt_line("Symbol* _morphology_sym;");
1323  }
1324 
1325  printer->pop_block(";");
1326 
1327  print_global_var_struct_assertions();
1328  print_global_var_struct_decl();
1329  print_global_var_external_access();
1330 
1331  print_global_param_default_values();
1332 }
1333 
1335  for (const auto& var: codegen_global_variables) {
1336  auto var_name = get_name(var);
1337  auto var_expr = get_variable_name(var_name, false);
1338 
1339  printer->fmt_push_block("auto {}() -> std::decay<decltype({})>::type ",
1340  method_name(var_name),
1341  var_expr);
1342  printer->fmt_line("return {};", var_expr);
1343  printer->pop_block();
1344  }
1345  if (!codegen_global_variables.empty()) {
1346  printer->add_newline();
1347  }
1348 
1349  for (const auto& var: info.external_variables) {
1350  auto var_name = get_name(var);
1351  printer->fmt_line("double {}();", var_name);
1352  }
1353  if (!info.external_variables.empty()) {
1354  printer->add_newline();
1355  }
1356 }
1357 
1359  printer->push_block("static std::vector<double> _parameter_defaults =");
1360 
1361  std::vector<std::string> defaults;
1362  for (const auto& p: info.range_parameter_vars) {
1363  double value = p->get_value() == nullptr ? 0.0 : *p->get_value();
1364  defaults.push_back(fmt::format("{:g} /* {} */", value, p->get_name()));
1365  }
1366 
1367  printer->add_multi_line(fmt::format("{}", fmt::join(defaults, ",\n")));
1368  printer->pop_block(";");
1369 }
1370 
1372  auto variable_printer = [&](const std::vector<SymbolType>& variables, bool if_array) {
1373  for (const auto& variable: variables) {
1374  if (variable->is_array() == if_array) {
1375  // false => do not use the instance struct, which is not
1376  // defined in the global declaration that we are printing
1377  auto name = get_variable_name(variable->get_name(), false);
1378  auto ename = add_escape_quote(variable->get_name() + "_" + info.mod_suffix);
1379  if (if_array) {
1380  auto length = variable->get_length();
1381  printer->fmt_line("{{{}, {}, {}}},", ename, name, length);
1382  } else {
1383  printer->fmt_line("{{{}, &{}}},", ename, name);
1384  }
1385  }
1386  }
1387  };
1388 
1389  auto globals = info.global_variables;
1390  auto thread_vars = info.thread_variables;
1391 
1392  if (info.table_count > 0) {
1393  globals.push_back(make_symbol(naming::USE_TABLE_VARIABLE));
1394  }
1395 
1396  printer->add_newline(2);
1397  printer->add_line("/** connect global (scalar) variables to hoc -- */");
1398  printer->add_line("static DoubScal hoc_scalar_double[] = {");
1399  printer->increase_indent();
1400  variable_printer(globals, false);
1401  variable_printer(thread_vars, false);
1402  printer->add_line("{nullptr, nullptr}");
1403  printer->decrease_indent();
1404  printer->add_line("};");
1405 
1406  printer->add_newline(2);
1407  printer->add_line("/** connect global (array) variables to hoc -- */");
1408  printer->add_line("static DoubVec hoc_vector_double[] = {");
1409  printer->increase_indent();
1410  variable_printer(globals, true);
1411  variable_printer(thread_vars, true);
1412  printer->add_line("{nullptr, nullptr, 0}");
1413  printer->decrease_indent();
1414  printer->add_line("};");
1415 
1416  printer->add_newline(2);
1417  printer->add_line("/* declaration of user functions */");
1418 
1419  auto print_entrypoint_decl = [this](const auto& callables, auto get_name) {
1420  for (const auto& node: callables) {
1421  const auto name = get_name(node);
1422  printer->fmt_line("{};", hoc_function_signature(name));
1423 
1424  if (!info.point_process) {
1425  printer->fmt_line("{};", py_function_signature(name));
1426  }
1427  }
1428  };
1429 
1430  auto get_name = [](const auto& node) { return node->get_node_name(); };
1431  print_entrypoint_decl(info.functions, get_name);
1432  print_entrypoint_decl(info.procedures, get_name);
1433  print_entrypoint_decl(info.function_tables, get_name);
1434  print_entrypoint_decl(info.function_tables, [](const auto& node) {
1435  auto node_name = node->get_node_name();
1436  return "table_" + node_name;
1437  });
1438 
1439  printer->add_newline(2);
1440  printer->add_line("/* connect user functions to hoc names */");
1441  printer->add_line("static VoidFunc hoc_intfunc[] = {");
1442  printer->increase_indent();
1443  if (info.point_process) {
1444  printer->add_line("{0, 0}");
1445  printer->decrease_indent();
1446  printer->add_line("};");
1447  printer->add_line("static Member_func _member_func[] = {");
1448  printer->increase_indent();
1449  printer->add_multi_line(R"CODE(
1450  {"loc", _hoc_loc_pnt},
1451  {"has_loc", _hoc_has_loc},
1452  {"get_loc", _hoc_get_loc_pnt},)CODE");
1453  } else {
1454  if (info.mod_suffix != "nothing") {
1455  printer->fmt_line("{{\"setdata_{}\", _hoc_setdata}},", info.mod_suffix);
1456  }
1457  }
1458 
1459  auto print_callable_reg = [this](const auto& callables, auto get_name) {
1460  for (const auto& node: callables) {
1461  const auto name = get_name(node);
1462  printer->fmt_line("{{\"{}{}\", {}}},", name, info.rsuffix, hoc_function_name(name));
1463  }
1464  };
1465 
1466  print_callable_reg(info.procedures, get_name);
1467  print_callable_reg(info.functions, get_name);
1468  print_callable_reg(info.function_tables, get_name);
1469  print_callable_reg(info.function_tables, [](const auto& node) {
1470  auto node_name = node->get_node_name();
1471  return "table_" + node_name;
1472  });
1473 
1474  printer->add_line("{nullptr, nullptr}");
1475  printer->decrease_indent();
1476  printer->add_line("};");
1477 
1478 
1479  auto print_py_callable_reg = [this](const auto& callables, auto get_name) {
1480  for (const auto& callable: callables) {
1481  const auto name = get_name(callable);
1482  printer->fmt_line("{{\"{}\", {}}},", name, py_function_name(name));
1483  }
1484  };
1485 
1486  if (!info.point_process) {
1487  printer->push_block("static NPyDirectMechFunc npy_direct_func_proc[] =");
1488  print_py_callable_reg(info.procedures,
1489  [](const auto& callable) { return callable->get_node_name(); });
1490  print_py_callable_reg(info.functions,
1491  [](const auto& callable) { return callable->get_node_name(); });
1492  print_py_callable_reg(info.function_tables,
1493  [](const auto& callable) { return callable->get_node_name(); });
1494  print_py_callable_reg(info.function_tables, [](const auto& callable) {
1495  return "table_" + callable->get_node_name();
1496  });
1497  printer->add_line("{nullptr, nullptr}");
1498  printer->pop_block(";");
1499  }
1500 }
1501 
1503  printer->add_newline(2);
1504  printer->fmt_push_block("extern \"C\" void _{}_reg()", info.mod_file);
1505  if (info.mod_suffix == "nothing") {
1506  print_mechanism_register_nothing();
1507  } else {
1508  print_mechanism_register_regular();
1509  }
1510  printer->pop_block();
1511 }
1512 
1514  printer->add_line("_initlists();");
1515  printer->add_newline();
1516 
1517  for (const auto& ion: info.ions) {
1518  double valence = ion.valence.value_or(-10000.0);
1519  printer->fmt_line("ion_reg(\"{}\", {});", ion.name, valence);
1520  }
1521  if (!info.ions.empty()) {
1522  printer->add_newline();
1523  }
1524 
1525  for (const auto& ion: info.ions) {
1526  printer->fmt_line("_{0}_sym = hoc_lookup(\"{0}_ion\");", ion.name);
1527  }
1528  if (!info.ions.empty()) {
1529  printer->add_newline();
1530  }
1531 
1532  const auto compute_functions_parameters =
1533  breakpoint_exist()
1534  ? fmt::format("{}, {}, {}",
1535  nrn_cur_required() ? method_name(naming::NRN_CUR_METHOD) : "nullptr",
1536  method_name(naming::NRN_JACOB_METHOD),
1537  nrn_state_required() ? method_name(naming::NRN_STATE_METHOD) : "nullptr")
1538  : "nullptr, nullptr, nullptr";
1539 
1540 
1541  const auto register_mech_args = fmt::format("{}, {}, {}, {}, {}, {}",
1542  get_channel_info_var_name(),
1543  method_name(naming::NRN_ALLOC_METHOD),
1544  compute_functions_parameters,
1545  method_name(naming::NRN_INIT_METHOD),
1546  info.first_pointer_var_index,
1547  1 + info.thread_data_index);
1548  if (info.point_process) {
1549  printer->fmt_line(
1550  "_pointtype = point_register_mech({}, _hoc_create_pnt, _hoc_destroy_pnt, "
1551  "_member_func);",
1552  register_mech_args);
1553 
1554  if (info.destructor_node) {
1555  printer->fmt_line("register_destructor({});",
1556  method_name(naming::NRN_DESTRUCTOR_METHOD));
1557  }
1558  } else {
1559  printer->fmt_line("register_mech({});", register_mech_args);
1560  }
1561 
1562 
1563  if (info.thread_callback_register) {
1564  printer->fmt_line("_extcall_thread.resize({});", info.thread_data_index + 1);
1565  printer->fmt_line("thread_mem_init(_extcall_thread.data());");
1566  printer->fmt_line("{} = 0;", get_variable_name("thread_data_in_use", false));
1567  }
1568 
1569 
1570  /// type related information
1571  printer->add_newline();
1572  printer->fmt_line("mech_type = nrn_get_mechtype({}[1]);", get_channel_info_var_name());
1573 
1574  printer->add_line("hoc_register_parm_default(mech_type, &_parameter_defaults);");
1575 
1576  // register the table-checking function
1577  if (info.table_count > 0) {
1578  printer->fmt_line("_nrn_thread_table_reg(mech_type, {});", table_thread_function_name());
1579  }
1580 
1581  printer->add_line("_nrn_mechanism_register_data_fields(mech_type,");
1582  printer->increase_indent();
1583 
1584  const auto codegen_float_variables_size = codegen_float_variables.size();
1585  std::vector<std::string> mech_register_args;
1586 
1587  for (int i = 0; i < codegen_float_variables_size; ++i) {
1588  const auto& float_var = codegen_float_variables[i];
1589  if (float_var->is_array()) {
1590  mech_register_args.push_back(
1591  fmt::format("_nrn_mechanism_field<double>{{\"{}\", {}}} /* {} */",
1592  float_var->get_name(),
1593  float_var->get_length(),
1594  i));
1595  } else {
1596  mech_register_args.push_back(fmt::format(
1597  "_nrn_mechanism_field<double>{{\"{}\"}} /* {} */", float_var->get_name(), i));
1598  }
1599  }
1600 
1601  const auto codegen_int_variables_size = codegen_int_variables.size();
1602  for (int i = 0; i < codegen_int_variables_size; ++i) {
1603  const auto& int_var = codegen_int_variables[i];
1604  const auto& name = int_var.symbol->get_name();
1605  if (i != info.semantics[i].index) {
1606  throw std::runtime_error("Broken logic.");
1607  }
1608  const auto& semantic = info.semantics[i].name;
1609 
1610  auto type = "double*";
1612  type = "Point_process*";
1613  } else if (name == naming::TQITEM_VARIABLE) {
1614  type = "void*";
1615  } else if (stringutils::starts_with(name, "style_") &&
1616  stringutils::starts_with(semantic, "#") &&
1617  stringutils::ends_with(semantic, "_ion")) {
1618  type = "int*";
1619  } else if (semantic == naming::FOR_NETCON_SEMANTIC) {
1620  type = "void*";
1621  }
1622 
1623  mech_register_args.push_back(
1624  fmt::format("_nrn_mechanism_field<{}>{{\"{}\", \"{}\"}} /* {} */",
1625  type,
1626  name,
1627  info.semantics[i].name,
1628  i));
1629  }
1630 
1631  if (info.emit_cvode) {
1632  mech_register_args.push_back(
1633  fmt::format("_nrn_mechanism_field<int>{{\"{}\", \"cvodeieq\"}} /* {} */",
1635  codegen_int_variables_size));
1636  }
1637 
1638  printer->add_multi_line(fmt::format("{}", fmt::join(mech_register_args, ",\n")));
1639 
1640  printer->decrease_indent();
1641  printer->add_line(");");
1642  printer->add_newline();
1643 
1644  printer->fmt_line("hoc_register_prop_size(mech_type, {}, {});",
1645  float_variables_size(),
1646  int_variables_size() + static_cast<int>(info.emit_cvode));
1647 
1648  for (int i = 0; i < codegen_int_variables_size; ++i) {
1649  if (i != info.semantics[i].index) {
1650  throw std::runtime_error("Broken logic.");
1651  }
1652 
1653  printer->fmt_line("hoc_register_dparam_semantics(mech_type, {}, \"{}\");",
1654  i,
1655  info.semantics[i].name);
1656  }
1657 
1658  if (!info.longitudinal_diffusion_info.empty()) {
1659  printer->fmt_line("hoc_register_ldifus1(_apply_diffusion_function);");
1660  }
1661 
1662 
1663  if (info.write_concentration) {
1664  printer->fmt_line("nrn_writes_conc(mech_type, 0);");
1665  }
1666 
1667  if (info.artificial_cell) {
1668  printer->fmt_line("add_nrn_artcell(mech_type, {});", info.tqitem_index);
1669  }
1670 
1671  if (info.net_event_used) {
1672  printer->add_line("add_nrn_has_net_event(mech_type);");
1673  }
1674 
1675  if (info.for_netcon_used) {
1676  auto dparam_it =
1677  std::find_if(info.semantics.begin(), info.semantics.end(), [](const IndexSemantics& a) {
1678  return a.name == naming::FOR_NETCON_SEMANTIC;
1679  });
1680  if (dparam_it == info.semantics.end()) {
1681  throw std::runtime_error("Couldn't find `fornetcon` variable.");
1682  }
1683 
1684  int dparam_index = dparam_it->index;
1685  printer->fmt_line("add_nrn_fornetcons(mech_type, {});", dparam_index);
1686  }
1687 
1688  printer->add_line("hoc_register_var(hoc_scalar_double, hoc_vector_double, hoc_intfunc);");
1689  if (!info.point_process) {
1690  printer->add_line("hoc_register_npy_direct(mech_type, npy_direct_func_proc);");
1691  }
1692  if (info.net_receive_node) {
1693  printer->fmt_line("pnt_receive[mech_type] = nrn_net_receive_{};", info.mod_suffix);
1694  printer->fmt_line("pnt_receive_size[mech_type] = {};", info.num_net_receive_parameters);
1695  }
1696 
1697  if (info.net_receive_initial_node) {
1698  printer->fmt_line("pnt_receive_init[mech_type] = net_init;");
1699  }
1700 
1701  if (info.thread_callback_register) {
1702  printer->add_line("_nrn_thread_reg(mech_type, 1, thread_mem_init);");
1703  printer->add_line("_nrn_thread_reg(mech_type, 0, thread_mem_cleanup);");
1704  }
1705 
1706  if (info.diam_used) {
1707  printer->fmt_line("{}._morphology_sym = hoc_lookup(\"morphology\");",
1708  global_struct_instance());
1709  }
1710 
1711  if (info.emit_cvode) {
1712  printer->fmt_line("hoc_register_dparam_semantics(mech_type, {}, \"cvodeieq\");",
1713  codegen_int_variables_size);
1714  printer->fmt_line("hoc_register_cvode(mech_type, {}, {}, {}, {});",
1715  method_name(naming::CVODE_COUNT_NAME),
1718  method_name(naming::CVODE_SETUP_STIFF_NAME));
1719  printer->fmt_line("hoc_register_tolerance(mech_type, _hoc_state_tol, &_atollist);");
1720  }
1721 }
1722 
1723 
1725  printer->add_line("hoc_register_var(hoc_scalar_double, hoc_vector_double, hoc_intfunc);");
1726 }
1727 
1728 
1730  if (!info.thread_callback_register) {
1731  return;
1732  }
1733 
1734  auto static_thread_data = get_variable_name("thread_data", false);
1735  auto inuse = get_variable_name("thread_data_in_use", false);
1736  auto thread_data_index = info.thread_var_thread_id;
1737  printer->push_block("static void thread_mem_init(Datum* _thread) ");
1738  printer->push_block(fmt::format("if({})", inuse));
1739  printer->fmt_line("_thread[{}] = {{neuron::container::do_not_search, new double[{}]{{}}}};",
1741  info.thread_var_data_size + info.top_local_thread_size);
1742  printer->pop_block();
1743  printer->push_block("else");
1744  printer->fmt_line("_thread[{}] = {{neuron::container::do_not_search, {}}};",
1746  static_thread_data);
1747  printer->fmt_line("{} = 1;", inuse);
1748  printer->pop_block();
1749  printer->pop_block();
1750 
1751  printer->push_block("static void thread_mem_cleanup(Datum* _thread) ");
1752  printer->fmt_line("double * _thread_data_ptr = _thread[{}].get<double*>();", thread_data_index);
1753  printer->push_block(fmt::format("if(_thread_data_ptr == {})", static_thread_data));
1754  printer->fmt_line("{} = 0;", inuse);
1755  printer->pop_block();
1756  printer->push_block("else");
1757  printer->add_line("delete[] _thread_data_ptr;");
1758  printer->pop_block();
1759  printer->pop_block();
1760 }
1761 
1762 
1764  for (auto const& [var, type]: info.neuron_global_variables) {
1765  auto const name = var->get_name();
1766  printer->fmt_line("extern {} {};", type, name);
1767  }
1768 }
1769 
1771  auto const value_initialize = print_initializers ? "{}" : "";
1772  printer->add_newline(2);
1773  printer->add_line("/** all mechanism instance variables and global variables */");
1774  printer->fmt_push_block("struct {} ", instance_struct());
1775 
1776  for (auto const& [var, type]: info.neuron_global_variables) {
1777  auto const name = var->get_name();
1778  printer->fmt_line("{}* {}{};",
1779  type,
1780  name,
1781  print_initializers ? fmt::format("{{&::{}}}", name) : std::string{});
1782  }
1783  for (auto& var: codegen_float_variables) {
1784  const auto& name = var->get_name();
1785  printer->fmt_line("double* {}{};", name, value_initialize);
1786  }
1787  for (auto& var: codegen_int_variables) {
1788  const auto& name = var.symbol->get_name();
1789  auto position = position_of_int_var(name);
1790 
1792  continue;
1793  } else if (var.is_index || var.is_integer) {
1794  // In NEURON we don't create caches for `int*`. Hence, do nothing.
1795  } else if (info.semantics[position].name == naming::POINTER_SEMANTIC) {
1796  // we don't need these either.
1797  } else {
1798  auto qualifier = var.is_constant ? "const " : "";
1799  auto type = var.is_vdata ? "void*" : default_float_data_type();
1800  printer->fmt_line("{}{}* const* {}{};", qualifier, type, name, value_initialize);
1801  }
1802  }
1803 
1804  printer->fmt_line("{}* {}{};",
1805  global_struct(),
1807  print_initializers ? fmt::format("{{&{}}}", global_struct_instance())
1808  : std::string{});
1809  printer->pop_block(";");
1810 }
1811 
1813  printer->add_newline(2);
1814  printer->fmt_push_block("static {} make_instance_{}(_nrn_mechanism_cache_range* _lmc)",
1815  instance_struct(),
1816  info.mod_suffix);
1817 
1818  printer->push_block("if(_lmc == nullptr)");
1819  printer->fmt_line("return {}();", instance_struct());
1820  printer->pop_block_nl(2);
1821 
1822  printer->fmt_push_block("return {}", instance_struct());
1823 
1824  std::vector<std::string> make_instance_args;
1825 
1826 
1827  for (auto const& [var, type]: info.neuron_global_variables) {
1828  auto const name = var->get_name();
1829  make_instance_args.push_back(fmt::format("&::{}", name));
1830  }
1831 
1832 
1833  const auto codegen_float_variables_size = codegen_float_variables.size();
1834  for (int i = 0; i < codegen_float_variables_size; ++i) {
1835  const auto& float_var = codegen_float_variables[i];
1836  if (float_var->is_array()) {
1837  make_instance_args.push_back(
1838  fmt::format("_lmc->template data_array_ptr<{}, {}>()", i, float_var->get_length()));
1839  } else {
1840  make_instance_args.push_back(fmt::format("_lmc->template fpfield_ptr<{}>()", i));
1841  }
1842  }
1843 
1844  const auto codegen_int_variables_size = codegen_int_variables.size();
1845  for (size_t i = 0; i < codegen_int_variables_size; ++i) {
1846  const auto& var = codegen_int_variables[i];
1847  auto name = var.symbol->get_name();
1848  auto sem = info.semantics[i].name;
1849  auto const variable = [&var, &sem, i]() -> std::string {
1850  if (var.is_index || var.is_integer) {
1851  return "";
1852  } else if (var.is_vdata) {
1853  return "";
1854  } else if (sem == naming::POINTER_SEMANTIC) {
1855  return "";
1856  } else {
1857  return fmt::format("_lmc->template dptr_field_ptr<{}>()", i);
1858  }
1859  }();
1860  if (variable != "") {
1861  make_instance_args.push_back(variable);
1862  }
1863  }
1864 
1865  printer->add_multi_line(fmt::format("{}", fmt::join(make_instance_args, ",\n")));
1866 
1867  printer->pop_block(";");
1868  printer->pop_block();
1869 }
1870 
1872  printer->add_newline(2);
1873  printer->fmt_push_block("struct {} ", node_data_struct());
1874 
1875  // Pointers to node variables
1876  printer->add_line("int const * nodeindices;");
1877  printer->add_line("double const * node_voltages;");
1878  printer->add_line("double * node_diagonal;");
1879  printer->add_line("double * node_rhs;");
1880  printer->add_line("int nodecount;");
1881 
1882  printer->pop_block(";");
1883 }
1884 
1886  printer->add_newline(2);
1887  printer->fmt_push_block("static {} make_node_data_{}(NrnThread& nt, Memb_list& _ml_arg)",
1888  node_data_struct(),
1889  info.mod_suffix);
1890 
1891  std::vector<std::string> make_node_data_args = {"_ml_arg.nodeindices",
1892  "nt.node_voltage_storage()",
1893  "nt.node_d_storage()",
1894  "nt.node_rhs_storage()",
1895  "_ml_arg.nodecount"};
1896 
1897  printer->fmt_push_block("return {}", node_data_struct());
1898  printer->add_multi_line(fmt::format("{}", fmt::join(make_node_data_args, ",\n")));
1899 
1900  printer->pop_block(";");
1901  printer->pop_block();
1902 
1903 
1904  printer->fmt_push_block("static {} make_node_data_{}(Prop * _prop)",
1905  node_data_struct(),
1906  info.mod_suffix);
1907 
1908  printer->push_block("if(!_prop)");
1909  printer->fmt_line("return {}();", node_data_struct());
1910  printer->pop_block_nl(2);
1911 
1912  printer->add_line("static std::vector<int> node_index{0};");
1913  printer->add_line("Node* _node = _nrn_mechanism_access_node(_prop);");
1914 
1915  make_node_data_args = {"node_index.data()",
1916  "&_nrn_mechanism_access_voltage(_node)",
1917  "&_nrn_mechanism_access_d(_node)",
1918  "&_nrn_mechanism_access_rhs(_node)",
1919  "1"};
1920 
1921  printer->fmt_push_block("return {}", node_data_struct());
1922  printer->add_multi_line(fmt::format("{}", fmt::join(make_node_data_args, ",\n")));
1923 
1924  printer->pop_block(";");
1925  printer->pop_block();
1926  printer->add_newline();
1927 }
1928 
1930  if (codegen_thread_variables.empty()) {
1931  return;
1932  }
1933 
1934  printer->add_newline(2);
1935  printer->fmt_push_block("struct {} ", thread_variables_struct());
1936  printer->add_line("double * thread_data;");
1937  printer->add_newline();
1938 
1939  std::string simd_width = "1";
1940 
1941 
1942  for (const auto& var_info: codegen_thread_variables) {
1943  printer->fmt_push_block("double * {}_ptr(size_t id)", var_info.symbol->get_name());
1944  printer->fmt_line("return thread_data + {} + (id % {});", var_info.offset, simd_width);
1945  printer->pop_block();
1946 
1947  printer->fmt_push_block("double & {}(size_t id)", var_info.symbol->get_name());
1948  printer->fmt_line("return thread_data[{} + (id % {})];", var_info.offset, simd_width);
1949  printer->pop_block();
1950  }
1951  printer->add_newline();
1952 
1953  printer->push_block(fmt::format("{}(double * const thread_data)", thread_variables_struct()));
1954  printer->fmt_line("this->thread_data = thread_data;");
1955  printer->pop_block();
1956 
1957  printer->pop_block(";");
1958 }
1959 
1960 
1962  // read ion statements
1963  auto read_statements = ion_read_statements(BlockType::Initial);
1964  for (auto& statement: read_statements) {
1965  printer->add_line(statement);
1966  }
1967 
1968  // initial block
1969  if (node != nullptr) {
1970  const auto& block = node->get_statement_block();
1971  print_statement_block(*block, false, false);
1972  }
1973 
1974  // write ion statements
1975  auto write_statements = ion_write_statements(BlockType::Initial);
1976  for (auto& statement: write_statements) {
1977  auto text = process_shadow_update_statement(statement, BlockType::Initial);
1978  printer->add_line(text);
1979  }
1980 }
1981 
1983  if (info.mod_suffix == "nothing") {
1984  return;
1985  }
1986 
1987  printer->add_line(
1988  "_nrn_mechanism_cache_range _lmc{_sorted_token, *nt, *_ml_arg, _ml_arg->type()};");
1989  printer->fmt_line("auto inst = make_instance_{}(&_lmc);", info.mod_suffix);
1990  if (!info.artificial_cell) {
1991  printer->fmt_line("auto node_data = make_node_data_{}(*nt, *_ml_arg);", info.mod_suffix);
1992  }
1993  printer->add_line("auto* _thread = _ml_arg->_thread;");
1994  if (!codegen_thread_variables.empty()) {
1995  printer->fmt_line("auto _thread_vars = {}(_thread[{}].get<double*>());",
1996  thread_variables_struct(),
1997  info.thread_var_thread_id);
1998  }
1999 }
2000 
2001 
2003  if (info.mod_suffix == "nothing") {
2004  return;
2005  }
2006 
2007  printer->add_line("Datum* _ppvar = _nrn_mechanism_access_dparam(prop);");
2008  printer->add_line("_nrn_mechanism_cache_instance _lmc{prop};");
2009  printer->add_line("const size_t id = 0;");
2010 
2011  printer->fmt_line("auto inst = make_instance_{}(prop ? &_lmc : nullptr);", info.mod_suffix);
2012  if (!info.artificial_cell) {
2013  printer->fmt_line("auto node_data = make_node_data_{}(prop);", info.mod_suffix);
2014  }
2015 
2016  if (!codegen_thread_variables.empty()) {
2017  printer->fmt_line("auto _thread_vars = {}({}_global.thread_data);",
2018  thread_variables_struct(),
2019  info.mod_suffix);
2020  }
2021 
2022  printer->add_newline();
2023 }
2024 
2025 
2027  const std::string& function_name) {
2028  std::string method = function_name.empty() ? compute_method_name(type) : function_name;
2029  ParamVector args = {{"", "const _nrn_model_sorted_token&", "", "_sorted_token"},
2030  {"", "NrnThread*", "", "nt"},
2031  {"", "Memb_list*", "", "_ml_arg"},
2032  {"", "int", "", "_type"}};
2033  printer->fmt_push_block("static void {}({})", method, get_parameter_str(args));
2034  print_entrypoint_setup_code_from_memb_list();
2035  printer->add_line("auto nodecount = _ml_arg->nodecount;");
2036 }
2037 
2038 
2039 void CodegenNeuronCppVisitor::print_nrn_init(bool skip_init_check) {
2040  printer->add_newline(2);
2041 
2042  print_global_function_common_code(BlockType::Initial);
2043 
2044  printer->push_block("for (int id = 0; id < nodecount; id++)");
2045 
2046  printer->add_line("auto* _ppvar = _ml_arg->pdata[id];");
2047  if (!info.artificial_cell) {
2048  printer->add_line("int node_id = node_data.nodeindices[id];");
2049  printer->add_line("inst.v_unused[id] = node_data.node_voltages[node_id];");
2050  }
2051 
2052  print_rename_state_vars();
2053 
2054  if (!info.changed_dt.empty()) {
2055  printer->fmt_line("double _save_prev_dt = {};",
2056  get_variable_name(naming::NTHREAD_DT_VARIABLE));
2057  printer->fmt_line("{} = {};",
2058  get_variable_name(naming::NTHREAD_DT_VARIABLE),
2059  info.changed_dt);
2060  }
2061 
2062  print_initial_block(info.initial_node);
2063 
2064  if (!info.changed_dt.empty()) {
2065  printer->fmt_line("{} = _save_prev_dt;", get_variable_name(naming::NTHREAD_DT_VARIABLE));
2066  }
2067 
2068  printer->pop_block();
2069  printer->pop_block();
2070 }
2071 
2073  printer->add_newline(2);
2074 
2075  ParamVector args = {{"", "const _nrn_model_sorted_token&", "", "_sorted_token"},
2076  {"", "NrnThread*", "", "nt"},
2077  {"", "Memb_list*", "", "_ml_arg"},
2078  {"", "int", "", "_type"}};
2079 
2080  printer->fmt_push_block("static void {}({})",
2081  method_name(naming::NRN_JACOB_METHOD),
2082  get_parameter_str(args)); // begin function
2083 
2084 
2085  print_entrypoint_setup_code_from_memb_list();
2086  printer->fmt_line("auto nodecount = _ml_arg->nodecount;");
2087  printer->push_block("for (int id = 0; id < nodecount; id++)"); // begin for
2088 
2089  if (breakpoint_exist()) {
2090  printer->add_line("int node_id = node_data.nodeindices[id];");
2091  printer->fmt_line("node_data.node_diagonal[node_id] {} inst.{}[id];",
2092  operator_for_d(),
2095  }
2096 
2097  printer->pop_block(); // end for
2098  printer->pop_block(); // end function
2099 }
2100 
2101 
2103  if (info.constructor_node) {
2104  printer->fmt_line("static void {}(Prop* prop);",
2105  method_name(naming::NRN_CONSTRUCTOR_METHOD));
2106  }
2107 }
2108 
2110  if (info.constructor_node) {
2111  printer->fmt_push_block("static void {}(Prop* prop)",
2112  method_name(naming::NRN_CONSTRUCTOR_METHOD));
2113 
2114  print_entrypoint_setup_code_from_prop();
2115 
2116  auto block = info.constructor_node->get_statement_block();
2117  print_statement_block(*block, false, false);
2118 
2119  printer->pop_block();
2120  }
2121 }
2122 
2123 
2125  printer->fmt_line("static void {}(Prop* prop);", method_name(naming::NRN_DESTRUCTOR_METHOD));
2126 }
2127 
2129  printer->fmt_push_block("static void {}(Prop* prop)",
2130  method_name(naming::NRN_DESTRUCTOR_METHOD));
2131  print_entrypoint_setup_code_from_prop();
2132 
2133  for (const auto& rv: info.random_variables) {
2134  printer->fmt_line("nrnran123_deletestream((nrnran123_State*) {});",
2135  get_variable_name(get_name(rv), false));
2136  }
2137 
2138 
2139  if (info.destructor_node) {
2140  auto block = info.destructor_node->get_statement_block();
2141  print_statement_block(*block, false, false);
2142  }
2143 
2144  printer->pop_block();
2145 }
2146 
2147 
2149  printer->add_newline(2);
2150 
2151  auto method = method_name(naming::NRN_ALLOC_METHOD);
2152  printer->fmt_push_block("static void {}(Prop* _prop)", method);
2153  printer->add_line("Datum *_ppvar = nullptr;");
2154 
2155  if (info.point_process) {
2156  printer->push_block("if (nrn_point_prop_)");
2157  printer->add_multi_line(R"CODE(
2158  _nrn_mechanism_access_alloc_seq(_prop) = _nrn_mechanism_access_alloc_seq(nrn_point_prop_);
2159  _ppvar = _nrn_mechanism_access_dparam(nrn_point_prop_);
2160  )CODE");
2161  printer->chain_block("else");
2162  }
2163  if (info.semantic_variable_count || info.emit_cvode) {
2164  printer->fmt_line("_ppvar = nrn_prop_datum_alloc(mech_type, {}, _prop);",
2165  info.semantic_variable_count + static_cast<int>(info.emit_cvode));
2166  printer->add_line("_nrn_mechanism_access_dparam(_prop) = _ppvar;");
2167  }
2168  printer->add_multi_line(R"CODE(
2169  _nrn_mechanism_cache_instance _lmc{_prop};
2170  size_t const _iml = 0;
2171  )CODE");
2172  printer->fmt_line("assert(_nrn_mechanism_get_num_vars(_prop) == {});",
2173  codegen_float_variables.size());
2174  if (float_variables_size()) {
2175  printer->add_line("/*initialize range parameters*/");
2176  for (size_t i_param = 0; i_param < info.range_parameter_vars.size(); ++i_param) {
2177  const auto var = info.range_parameter_vars[i_param];
2178  if (var->is_array()) {
2179  continue;
2180  }
2181  const auto& var_name = var->get_name();
2182  auto var_pos = position_of_float_var(var_name);
2183 
2184  printer->fmt_line("_lmc.template fpfield<{}>(_iml) = {}; /* {} */",
2185  var_pos,
2186  fmt::format("_parameter_defaults[{}]", i_param),
2187  var_name);
2188  }
2189  }
2190  if (info.point_process) {
2191  printer->pop_block();
2192  }
2193 
2194  if (info.semantic_variable_count) {
2195  printer->add_line("_nrn_mechanism_access_dparam(_prop) = _ppvar;");
2196  }
2197 
2198  const auto codegen_int_variables_size = codegen_int_variables.size();
2199 
2200  if (info.diam_used || info.area_used) {
2201  for (size_t i = 0; i < codegen_int_variables.size(); ++i) {
2202  auto var_info = codegen_int_variables[i];
2203  if (var_info.symbol->get_name() == naming::DIAM_VARIABLE) {
2204  printer->fmt_line("Prop * morphology_prop = need_memb({}._morphology_sym);",
2205  global_struct_instance());
2206  printer->fmt_line(
2207  "_ppvar[{}] = _nrn_mechanism_get_param_handle(morphology_prop, 0);", i);
2208  }
2209  if (var_info.symbol->get_name() == naming::AREA_VARIABLE) {
2210  printer->fmt_line("_ppvar[{}] = _nrn_mechanism_get_area_handle(nrn_alloc_node_);",
2211  i);
2212  }
2213  }
2214  }
2215 
2216  for (const auto& ion: info.ions) {
2217  printer->fmt_line("Symbol * {}_sym = hoc_lookup(\"{}_ion\");", ion.name, ion.name);
2218  printer->fmt_line("Prop * {}_prop = need_memb({}_sym);", ion.name, ion.name);
2219 
2220  if (ion.is_exterior_conc_written()) {
2221  printer->fmt_line("nrn_check_conc_write(_prop, {}_prop, 0);", ion.name);
2222  }
2223 
2224  if (ion.is_interior_conc_written()) {
2225  printer->fmt_line("nrn_check_conc_write(_prop, {}_prop, 1);", ion.name);
2226  }
2227 
2228  int conc = ion.is_conc_written() ? 3 : int(ion.is_conc_read());
2229  int rev = ion.is_rev_written() ? 3 : int(ion.is_rev_read());
2230 
2231  printer->fmt_line("nrn_promote({}_prop, {}, {});", ion.name, conc, rev);
2232 
2233  for (size_t i = 0; i < codegen_int_variables_size; ++i) {
2234  const auto& var = codegen_int_variables[i];
2235 
2236  const std::string& var_name = var.symbol->get_name();
2237 
2238  if (stringutils::starts_with(var_name, "ion_")) {
2239  std::string ion_var_name = std::string(var_name.begin() + 4, var_name.end());
2240  if (ion.is_ionic_variable(ion_var_name) ||
2241  ion.is_current_derivative(ion_var_name) || ion.is_rev_potential(ion_var_name)) {
2242  printer->fmt_line("_ppvar[{}] = _nrn_mechanism_get_param_handle({}_prop, {});",
2243  i,
2244  ion.name,
2245  ion.variable_index(ion_var_name));
2246  }
2247  } else {
2248  if (ion.is_style(var_name)) {
2249  printer->fmt_line(
2250  "_ppvar[{}] = {{neuron::container::do_not_search, "
2251  "&(_nrn_mechanism_access_dparam({}_prop)[0].literal_value<int>())}};",
2252  i,
2253  ion.name);
2254  }
2255  }
2256  }
2257  }
2258 
2259  if (!info.random_variables.empty()) {
2260  for (const auto& rv: info.random_variables) {
2261  auto position = position_of_int_var(get_name(rv));
2262  printer->fmt_line("_ppvar[{}].literal_value<void*>() = nrnran123_newstream();",
2263  position);
2264  }
2265  printer->fmt_line("nrn_mech_inst_destruct[mech_type] = neuron::{};",
2266  method_name(naming::NRN_DESTRUCTOR_METHOD));
2267  }
2268 
2269  if (info.point_process || info.artificial_cell) {
2270  printer->fmt_push_block("if(!nrn_point_prop_)");
2271 
2272  if (info.constructor_node) {
2273  printer->fmt_line("{}(_prop);", method_name(naming::NRN_CONSTRUCTOR_METHOD));
2274  }
2275  printer->pop_block();
2276  }
2277 
2278  printer->pop_block();
2279 }
2280 
2281 
2282 /****************************************************************************************/
2283 /* Print nrn_state routine */
2284 /****************************************************************************************/
2285 
2287  if (!nrn_state_required()) {
2288  return;
2289  }
2290 
2291  printer->add_newline(2);
2292  print_global_function_common_code(BlockType::State);
2293 
2294  printer->push_block("for (int id = 0; id < nodecount; id++)");
2295  printer->add_line("int node_id = node_data.nodeindices[id];");
2296  printer->add_line("auto* _ppvar = _ml_arg->pdata[id];");
2297  if (!info.artificial_cell) {
2298  printer->add_line("inst.v_unused[id] = node_data.node_voltages[node_id];");
2299  }
2300 
2301  /**
2302  * \todo Eigen solver node also emits IonCurVar variable in the functor
2303  * but that shouldn't update ions in derivative block
2304  */
2305  if (ion_variable_struct_required()) {
2306  throw std::runtime_error("Not implemented.");
2307  }
2308 
2309  auto read_statements = ion_read_statements(BlockType::State);
2310  for (auto& statement: read_statements) {
2311  printer->add_line(statement);
2312  }
2313 
2314  if (info.nrn_state_block) {
2315  info.nrn_state_block->visit_children(*this);
2316  }
2317 
2318  for (auto& block: info.matexp_blocks) {
2319  if (!block->get_steadystate().get()->eval()) {
2320  block->accept(*this);
2321  }
2322  }
2323 
2324  if (info.currents.empty() && info.breakpoint_node != nullptr) {
2325  auto block = info.breakpoint_node->get_statement_block();
2326  print_statement_block(*block, false, false);
2327  }
2328 
2329  const auto& write_statements = ion_write_statements(BlockType::State);
2330  for (auto& statement: write_statements) {
2331  const auto& text = process_shadow_update_statement(statement, BlockType::State);
2332  printer->add_line(text);
2333  }
2334 
2335  printer->pop_block();
2336  printer->pop_block();
2337 }
2338 
2339 
2340 /****************************************************************************************/
2341 /* Print nrn_cur related routines */
2342 /****************************************************************************************/
2343 
2345  return get_arg_str(nrn_current_parameters());
2346 }
2347 
2348 
2350  if (ion_variable_struct_required()) {
2351  throw std::runtime_error("Not implemented.");
2352  }
2353 
2354  ParamVector params = {{"", "_nrn_mechanism_cache_range&", "", "_lmc"},
2355  {"", "NrnThread*", "", "nt"},
2356  {"", "Datum*", "", "_ppvar"},
2357  {"", "Datum*", "", "_thread"}};
2358 
2359  if (info.thread_callback_register) {
2360  auto type_name = fmt::format("{}&", thread_variables_struct());
2361  params.emplace_back("", type_name, "", "_thread_vars");
2362  }
2363  params.emplace_back("", "size_t", "", "id");
2364  params.emplace_back("", fmt::format("{}&", instance_struct()), "", "inst");
2365  params.emplace_back("", fmt::format("{}&", node_data_struct()), "", "node_data");
2366  params.emplace_back("", "double", "", "v");
2367  return params;
2368 }
2369 
2370 
2372  const auto& args = nrn_current_parameters();
2373  const auto& block = node.get_statement_block();
2374  printer->add_newline(2);
2375  printer->fmt_push_block("static inline double nrn_current_{}({})",
2376  info.mod_suffix,
2377  get_parameter_str(args));
2378  printer->add_line("inst.v_unused[id] = v;");
2379  printer->add_line("double current = 0.0;");
2380  print_statement_block(*block, false, false);
2381  for (auto& current: info.currents) {
2382  const auto& name = get_variable_name(current);
2383  printer->fmt_line("current += {};", name);
2384  }
2385  printer->add_line("return current;");
2386  printer->pop_block();
2387 }
2388 
2389 
2391  const auto& block = node.get_statement_block();
2392  print_statement_block(*block, false, false);
2393  if (!info.currents.empty()) {
2394  std::string sum;
2395  for (const auto& current: info.currents) {
2396  auto var = breakpoint_current(current);
2397  sum += get_variable_name(var);
2398  if (&current != &info.currents.back()) {
2399  sum += "+";
2400  }
2401  }
2402  printer->fmt_line("double rhs = {};", sum);
2403  }
2404 
2405  std::string sum;
2406  for (const auto& conductance: info.conductances) {
2407  auto var = breakpoint_current(conductance.variable);
2408  sum += get_variable_name(var);
2409  if (&conductance != &info.conductances.back()) {
2410  sum += "+";
2411  }
2412  }
2413  printer->fmt_line("double g = {};", sum);
2414 
2415  for (const auto& conductance: info.conductances) {
2416  if (!conductance.ion.empty()) {
2417  const auto& lhs = std::string(naming::ION_VARNAME_PREFIX) + "di" + conductance.ion +
2418  "dv";
2419  const auto& rhs = get_variable_name(conductance.variable);
2420  const ShadowUseStatement statement{lhs, "+=", rhs};
2421  const auto& text = process_shadow_update_statement(statement, BlockType::Equation);
2422  printer->add_line(text);
2423  }
2424  }
2425 }
2426 
2427 
2429  printer->fmt_line("double I1 = nrn_current_{}({}+0.001);",
2430  info.mod_suffix,
2431  nrn_current_arguments());
2432  for (auto& ion: info.ions) {
2433  for (auto& var: ion.writes) {
2434  if (ion.is_ionic_current(var)) {
2435  const auto& name = get_variable_name(var);
2436  printer->fmt_line("double di{} = {};", ion.name, name);
2437  }
2438  }
2439  }
2440  printer->fmt_line("double I0 = nrn_current_{}({});", info.mod_suffix, nrn_current_arguments());
2441  printer->add_line("double rhs = I0;");
2442 
2443  printer->add_line("double g = (I1-I0)/0.001;");
2444  for (auto& ion: info.ions) {
2445  for (auto& var: ion.writes) {
2446  if (ion.is_ionic_current(var)) {
2447  const auto& lhs = std::string(naming::ION_VARNAME_PREFIX) + "di" + ion.name + "dv";
2448  auto rhs = fmt::format("(di{}-{})/0.001", ion.name, get_variable_name(var));
2449  if (info.point_process) {
2450  auto area = get_variable_name(naming::NODE_AREA_VARIABLE);
2451  rhs += fmt::format("*1.e2/{}", area);
2452  }
2453  const ShadowUseStatement statement{lhs, "+=", rhs};
2454  const auto& text = process_shadow_update_statement(statement, BlockType::Equation);
2455  printer->add_line(text);
2456  }
2457  }
2458  }
2459 }
2460 
2461 
2463  printer->add_line("int node_id = node_data.nodeindices[id];");
2464  printer->add_line("double v = node_data.node_voltages[node_id];");
2465  printer->add_line("auto* _ppvar = _ml_arg->pdata[id];");
2466  const auto& read_statements = ion_read_statements(BlockType::Equation);
2467  for (auto& statement: read_statements) {
2468  printer->add_line(statement);
2469  }
2470 
2471  if (info.conductances.empty()) {
2472  print_nrn_cur_non_conductance_kernel();
2473  } else {
2474  print_nrn_cur_conductance_kernel(node);
2475  }
2476 
2477  const auto& write_statements = ion_write_statements(BlockType::Equation);
2478  for (auto& statement: write_statements) {
2479  auto text = process_shadow_update_statement(statement, BlockType::Equation);
2480  printer->add_line(text);
2481  }
2482 
2483  if (info.point_process) {
2484  const auto& area = get_variable_name(naming::NODE_AREA_VARIABLE);
2485  printer->fmt_line("double mfactor = 1.e2/{};", area);
2486  printer->add_line("g = g*mfactor;");
2487  printer->add_line("rhs = rhs*mfactor;");
2488  }
2489 
2490  // print_g_unused();
2491 }
2492 
2493 
2494 /// TODO: Edit for NEURON
2496  return;
2497 }
2498 
2499 
2500 /// TODO: Edit for NEURON
2502  if (!nrn_cur_required()) {
2503  return;
2504  }
2505 
2506  if (info.conductances.empty()) {
2507  print_nrn_current(*info.breakpoint_node);
2508  }
2509 
2510  printer->add_newline(2);
2511  printer->add_line("/** update current */");
2512  print_global_function_common_code(BlockType::Equation);
2513  // print_channel_iteration_block_parallel_hint(BlockType::Equation, info.breakpoint_node);
2514  printer->push_block("for (int id = 0; id < nodecount; id++)");
2515  print_nrn_cur_kernel(*info.breakpoint_node);
2516  // print_nrn_cur_matrix_shadow_update();
2517  // if (!nrn_cur_reduction_loop_required()) {
2518  // print_fast_imem_calculation();
2519  // }
2520 
2521 
2522  printer->fmt_line("node_data.node_rhs[node_id] {} rhs;", operator_for_rhs());
2523 
2524  if (breakpoint_exist()) {
2525  printer->fmt_line("inst.{}[id] = g;",
2528  }
2529  printer->pop_block();
2530 
2531  // if (nrn_cur_reduction_loop_required()) {
2532  // printer->push_block("for (int id = 0; id < nodecount; id++)");
2533  // print_nrn_cur_matrix_shadow_reduction();
2534  // printer->pop_block();
2535  // print_fast_imem_calculation();
2536  // }
2537 
2538  // print_kernel_data_present_annotation_block_end();
2539  printer->pop_block();
2540 }
2541 
2542 
2543 /****************************************************************************************/
2544 /* Main code printing entry points */
2545 /****************************************************************************************/
2546 
2548  print_standard_includes();
2549  print_neuron_includes();
2550 
2551  if (info.thread_callback_register) {
2552  printer->add_line("extern void _nrn_thread_reg(int, int, void(*)(Datum*));");
2553  }
2554 }
2555 
2556 
2558  print_global_macros();
2559  print_mechanism_variables_macros();
2560 
2561  printer->add_line("extern Node* nrn_alloc_node_;");
2562 }
2563 
2564 
2566  printer->add_newline();
2567  printer->add_line("/* NEURON global macro definitions */");
2568  if (info.vectorize) {
2569  printer->add_multi_line(R"CODE(
2570  /* VECTORIZED */
2571  #define NRN_VECTORIZED 1
2572  )CODE");
2573  } else {
2574  printer->add_multi_line(R"CODE(
2575  /* NOT VECTORIZED */
2576  #define NRN_VECTORIZED 0
2577  )CODE");
2578  }
2579 }
2580 
2581 
2583  printer->add_newline();
2584  printer->add_line("static constexpr auto number_of_datum_variables = ",
2585  std::to_string(int_variables_size()),
2586  ";");
2587  printer->add_line("static constexpr auto number_of_floating_point_variables = ",
2588  std::to_string(codegen_float_variables.size()),
2589  ";");
2590  printer->add_newline();
2591  printer->add_multi_line(R"CODE(
2592  namespace {
2593  template <typename T>
2594  using _nrn_mechanism_std_vector = std::vector<T>;
2595  using _nrn_model_sorted_token = neuron::model_sorted_token;
2596  using _nrn_mechanism_cache_range = neuron::cache::MechanismRange<number_of_floating_point_variables, number_of_datum_variables>;
2597  using _nrn_mechanism_cache_instance = neuron::cache::MechanismInstance<number_of_floating_point_variables, number_of_datum_variables>;
2598  using _nrn_non_owning_id_without_container = neuron::container::non_owning_identifier_without_container;
2599  template <typename T>
2600  using _nrn_mechanism_field = neuron::mechanism::field<T>;
2601  template <typename... Args>
2602  void _nrn_mechanism_register_data_fields(Args&&... args) {
2603  neuron::mechanism::register_data_fields(std::forward<Args>(args)...);
2604  }
2605  } // namespace
2606  )CODE");
2607 
2608  if (info.point_process) {
2609  printer->add_line("extern Prop* nrn_point_prop_;");
2610  } else {
2611  printer->add_line("Prop* hoc_getdata_range(int type);");
2612  }
2613  // for registration of tables
2614  if (info.table_count > 0) {
2615  printer->add_line("void _nrn_thread_table_reg(int, nrn_thread_table_check_t);");
2616  }
2617  // for CVODE
2618  printer->add_line("extern void _cvode_abstol(Symbol**, double*, int);");
2619  if (info.for_netcon_used) {
2620  printer->add_line("int _nrn_netcon_args(void*, double***);");
2621  }
2622 }
2623 
2624 
2625 void CodegenNeuronCppVisitor::print_data_structures(bool print_initializers) {
2626  print_mechanism_global_var_structure(print_initializers);
2627  print_mechanism_range_var_structure(print_initializers);
2628  print_node_data_structure(print_initializers);
2629  print_thread_variables_structure(print_initializers);
2630  print_make_instance();
2631  print_make_node_data();
2632 }
2633 
2634 
2636  if (!info.vectorize) {
2637  return;
2638  }
2639  printer->add_multi_line(R"CODE(
2640  #if NRN_PRCELLSTATE
2641  inst->v_unused[id] = v;
2642  #endif
2643  )CODE");
2644 }
2645 
2646 
2648  printer->add_multi_line(R"CODE(
2649  #if NRN_PRCELLSTATE
2650  inst->g_unused[id] = g;
2651  #endif
2652  )CODE");
2653 }
2654 
2656  print_hoc_py_wrapper_function_definitions();
2657  for (const auto& procedure: info.procedures) {
2658  print_procedure(*procedure);
2659  }
2660  for (const auto& function: info.functions) {
2661  print_function(*function);
2662  }
2663  for (const auto& function_table: info.function_tables) {
2664  print_function_tables(*function_table);
2665  }
2666 }
2667 
2669  print_nrn_init();
2670  print_nrn_cur();
2671  print_nrn_state();
2672  print_nrn_jacob();
2673  print_net_receive();
2674  print_net_init();
2675 }
2676 
2678  print_backend_info();
2679  print_headers_include();
2680  print_macro_definitions();
2681  print_neuron_global_variable_declarations();
2682  print_namespace_start();
2683  print_nmodl_constants();
2684  print_prcellstate_macros();
2685  print_mechanism_info();
2686  print_data_structures(true);
2687  print_nrn_constructor_declaration();
2688  print_nrn_destructor_declaration();
2689  print_nrn_alloc();
2690  print_function_prototypes();
2691  print_longitudinal_diffusion_callbacks();
2692  print_cvode_definitions();
2693  print_point_process_function_definitions();
2694  print_setdata_functions();
2695  print_check_table_entrypoint();
2696  print_top_verbatim_blocks();
2697  print_functors_definitions();
2698  print_global_variables_for_hoc();
2699  print_thread_memory_callbacks();
2700  print_function_definitions();
2701  print_compute_functions(); // only nrn_cur and nrn_state
2702  print_nrn_constructor();
2703  print_nrn_destructor();
2704  print_sdlists_init(true);
2705  print_mechanism_register();
2706  print_namespace_stop();
2707 }
2708 
2710  print_backend_info();
2711  print_headers_include();
2712  print_namespace_start();
2713  print_function_prototypes();
2714  print_top_verbatim_blocks();
2715  print_global_variables_for_hoc();
2716  print_function_definitions();
2717  print_mechanism_register();
2718  print_namespace_stop();
2719 }
2720 
2722  if (info.mod_suffix == "nothing") {
2723  print_codegen_routines_nothing();
2724  } else {
2725  print_codegen_routines_regular();
2726  }
2727 }
2728 
2730  throw std::runtime_error("Not implemented.");
2731 }
2732 
2733 
2735  auto const& arguments = node.get_arguments();
2736 
2737  if (printing_net_init) {
2738  throw std::runtime_error("Not implemented. [jfiwoei]");
2739  }
2740 
2741  std::string weight_pointer = "nullptr";
2742  auto point_process = get_variable_name(naming::POINT_PROCESS_VARIABLE,
2743  /* use_instance */ false);
2744  if (!printing_net_receive) {
2745  point_process += ".get<Point_process*>()";
2746  }
2747  const auto& tqitem = get_variable_name("tqitem", /* use_instance */ false);
2748 
2749  printer->fmt_text("{}(/* tqitem */ &{}, {}, {}, {} + ",
2750  info.artificial_cell ? "artcell_net_send" : "net_send",
2751  tqitem,
2752  weight_pointer,
2753  point_process,
2754  get_variable_name("t"));
2755  print_vector_elements(arguments, ", ");
2756  printer->add_text(')');
2757 }
2758 
2760  const auto& point_process = get_variable_name("point_process", /* use_instance */ false);
2761  const auto& tqitem = get_variable_name("tqitem", /* use_instance */ false);
2762 
2763  printer->fmt_text("{}(/* tqitem */ &{}, {}, ",
2764  info.artificial_cell ? "artcell_net_move" : "net_move",
2765  tqitem,
2766  point_process);
2767 
2768  print_vector_elements(node.get_arguments(), ", ");
2769  printer->add_text(')');
2770 }
2771 
2773  const auto& point_process = get_variable_name(naming::POINT_PROCESS_VARIABLE,
2774  /* use_instance */ false);
2775  printer->fmt_text("net_event({}, t)", point_process);
2776 }
2777 
2779  auto name = node.get_node_name();
2780  const auto& arguments = node.get_arguments();
2781  printer->add_text(method_name(name), "(");
2782  print_vector_elements(arguments, ", ");
2783  printer->add_text(')');
2784 }
2785 
2786 /**
2787  * Rename arguments to NET_RECEIVE block with corresponding pointer variable
2788  *
2789  * \code{.mod}
2790  * NET_RECEIVE (weight, R){
2791  * x = R
2792  * }
2793  * \endcode
2794  *
2795  * then generated code should be:
2796  *
2797  * \code{.cpp}
2798  * x[id] = _args[1];
2799  * \endcode
2800  *
2801  * So, the `R` in AST needs to be renamed with `_args[1]`.
2802  */
2803 static void rename_net_receive_arguments(const ast::NetReceiveBlock& net_receive_node,
2804  const ast::Node& node) {
2805  const auto& parameters = net_receive_node.get_parameters();
2806 
2807  auto n_parameters = parameters.size();
2808  for (size_t i = 0; i < n_parameters; ++i) {
2809  const auto& name = parameters[i]->get_node_name();
2810  auto var_used = VarUsageVisitor().variable_used(node, name);
2811  if (var_used) {
2812  RenameVisitor vr(name, fmt::format("_args[{}]", i));
2813  node.get_statement_block()->visit_children(vr);
2814  }
2815  }
2816 }
2817 
2818 
2820  return {{"", "Point_process*", "", "_pnt"},
2821  {"", "double*", "", "_args"},
2822  {"", "double", "", "flag"}};
2823 }
2824 
2825 
2827  printer->add_line("_nrn_mechanism_cache_instance _lmc{_pnt->prop};");
2828  printer->add_line("auto * nt = static_cast<NrnThread*>(_pnt->_vnt);");
2829  printer->add_line("auto * _ppvar = _nrn_mechanism_access_dparam(_pnt->prop);");
2830 
2831  printer->fmt_line("auto inst = make_instance_{}(&_lmc);", info.mod_suffix);
2832  if (!info.artificial_cell) {
2833  printer->fmt_line("auto node_data = make_node_data_{}(_pnt->prop);", info.mod_suffix);
2834  }
2835  printer->fmt_line("// nocmodl has a nullptr dereference for thread variables.");
2836  printer->fmt_line("// NMODL will fail to compile at a later point, because of");
2837  printer->fmt_line("// missing '_thread_vars'.");
2838  printer->fmt_line("Datum * _thread = nullptr;");
2839 
2840  printer->add_line("size_t id = 0;");
2841  printer->add_line("double t = nt->_t;");
2842 }
2843 
2845  return {{"", "const _nrn_model_sorted_token&", "", "_sorted_token"},
2846  {"", "NrnThread*", "", "nt"},
2847  {"", "Memb_list*", "", "_ml_arg"},
2848  {"", "int", "", "_type"}};
2849 }
2850 
2852  ParamVector args = {{"", "_nrn_mechanism_cache_range&", "", "_lmc"},
2853  {"", fmt::format("{}_Instance&", info.mod_suffix), "", "inst"},
2854  {"", fmt::format("{}_NodeData&", info.mod_suffix), "", "node_data"},
2855  {"", "size_t", "", "id"},
2856  {"", "Datum*", "", "_ppvar"},
2857  {"", "Datum*", "", "_thread"},
2858  {"", "NrnThread*", "", "nt"}};
2859 
2860  if (info.thread_callback_register) {
2861  auto type_name = fmt::format("{}&", thread_variables_struct());
2862  args.emplace_back("", type_name, "", "_thread_vars");
2863  }
2864  return args;
2865 }
2866 
2867 
2869  printer->fmt_push_block("static constexpr int {}(int _type)",
2870  method_name(naming::CVODE_COUNT_NAME));
2871  printer->fmt_line("return {};", info.cvode_block->get_n_odes()->get_value());
2872  printer->pop_block();
2873  printer->add_newline(2);
2874 }
2875 
2877  const CodegenNeuronCppVisitor::ParamVector tolerances_parameters = {
2878  {"", "Prop*", "", "_prop"},
2879  {"", "int", "", "equation_index"},
2880  {"", "neuron::container::data_handle<double>*", "", "_pv"},
2881  {"", "neuron::container::data_handle<double>*", "", "_pvdot"},
2882  {"", "double*", "", "_atol"},
2883  {"", "int", "", "_type"}};
2884 
2885  auto get_param_name = [](const auto& item) { return std::get<3>(item); };
2886 
2887  auto prop_name = get_param_name(tolerances_parameters[0]);
2888  auto eqindex_name = get_param_name(tolerances_parameters[1]);
2889  auto pv_name = get_param_name(tolerances_parameters[2]);
2890  auto pvdot_name = get_param_name(tolerances_parameters[3]);
2891  auto atol_name = get_param_name(tolerances_parameters[4]);
2892 
2893  printer->fmt_push_block("static void {}({})",
2895  get_parameter_str(tolerances_parameters));
2896  printer->fmt_line("auto* _ppvar = _nrn_mechanism_access_dparam({});", prop_name);
2897  printer->fmt_line("_ppvar[{}].literal_value<int>() = {};", int_variables_size(), eqindex_name);
2898  printer->fmt_push_block("for (int i = 0; i < {}(0); i++)",
2899  method_name(naming::CVODE_COUNT_NAME));
2900  printer->fmt_line("{}[i] = _nrn_mechanism_get_param_handle({}, _slist1[i]);",
2901  pv_name,
2902  prop_name);
2903  printer->fmt_line("{}[i] = _nrn_mechanism_get_param_handle({}, _dlist1[i]);",
2904  pvdot_name,
2905  prop_name);
2906  printer->fmt_line("_cvode_abstol(_atollist, {}, i);", atol_name);
2907  printer->pop_block();
2908  printer->pop_block();
2909  printer->add_newline(2);
2910 }
2911 
2913  const ast::StatementBlock& block) {
2914  printer->fmt_push_block("static int {}({})",
2915  method_name(name),
2916  get_parameter_str(cvode_update_parameters()));
2917  printer->add_line(
2918  "auto v = node_data.node_voltages ? "
2919  "node_data.node_voltages[node_data.nodeindices[id]] : 0.0;");
2920 
2921  print_statement_block(block, false, false);
2922 
2923  printer->add_line("return 0;");
2924  printer->pop_block();
2925  printer->add_newline(2);
2926 }
2927 
2928 void CodegenNeuronCppVisitor::print_cvode_setup(const std::string& setup_name,
2929  const std::string& update_name) {
2930  printer->fmt_push_block("static void {}({})",
2931  method_name(setup_name),
2932  get_parameter_str(cvode_setup_parameters()));
2933  print_entrypoint_setup_code_from_memb_list();
2934  printer->fmt_line("auto nodecount = _ml_arg->nodecount;");
2935  printer->push_block("for (int id = 0; id < nodecount; id++)");
2936  printer->add_line("auto* _ppvar = _ml_arg->pdata[id];");
2937  printer->add_line(
2938  "auto v = node_data.node_voltages ? "
2939  "node_data.node_voltages[node_data.nodeindices[id]] : 0.0;");
2940  printer->fmt_line("{}({});", method_name(update_name), get_arg_str(cvode_update_parameters()));
2941 
2942  printer->pop_block();
2943  printer->pop_block();
2944 
2945  printer->add_newline(2);
2946 }
2947 
2949  if (!info.emit_cvode) {
2950  return;
2951  }
2952 
2953  printer->add_newline(2);
2954  printer->add_line("/* Functions related to CVODE codegen */");
2955  print_cvode_count();
2956  print_cvode_tolerances();
2957  print_cvode_update(naming::CVODE_UPDATE_NON_STIFF_NAME,
2958  *info.cvode_block->get_non_stiff_block());
2959  print_cvode_update(naming::CVODE_UPDATE_STIFF_NAME, *info.cvode_block->get_stiff_block());
2962 }
2963 
2965  printing_net_receive = true;
2966  auto node = info.net_receive_node;
2967  if (!node) {
2968  return;
2969  }
2970 
2971  printer->fmt_push_block("static void nrn_net_receive_{}({})",
2972  info.mod_suffix,
2973  get_parameter_str(net_receive_args()));
2974 
2976  print_net_receive_common_code();
2977 
2978 
2979  print_statement_block(*node->get_statement_block(), false, false);
2980 
2981  printer->add_newline();
2982  printer->pop_block();
2983  printing_net_receive = false;
2984 }
2985 
2987  const auto node = info.net_receive_initial_node;
2988  if (node == nullptr) {
2989  return;
2990  }
2991 
2992  // rename net_receive arguments used in the initial block of net_receive
2993  rename_net_receive_arguments(*info.net_receive_node, *node);
2994 
2995  printing_net_init = true;
2996  printer->add_newline(2);
2997  printer->fmt_push_block("static void net_init({})", get_parameter_str(net_receive_args()));
2998 
2999  auto block = node->get_statement_block().get();
3000  if (!block->get_statements().empty()) {
3001  print_net_receive_common_code();
3002  print_statement_block(*block, false, false);
3003  }
3004  printer->pop_block();
3005  printing_net_init = false;
3006 }
3007 
3008 
3009 /****************************************************************************************/
3010 /* Overloaded visitor routines */
3011 /****************************************************************************************/
3012 
3013 /// TODO: Edit for NEURON
3015  return;
3016 }
3017 
3019  const ast::LongitudinalDiffusionBlock& /* node */) {
3020  // These are handled via `print_longitdudinal_*`.
3021 }
3022 
3024  // These are handled via `print_longitdudinal_*`.
3025 }
3026 
3028  // The setup for enabling this loop is:
3029  // double ** _fornetcon_data = ...;
3030  // for(size_t i = 0; i < n_netcons; ++i) {
3031  // double * _netcon_data = _fornetcon_data[i];
3032  //
3033  // // loop body.
3034  // }
3035  //
3036  // Where `_fornetcon_data` is an array of pointers to the arguments, one
3037  // for each netcon; and `_netcon_data` points to the arguments for the
3038  // current netcon.
3039  //
3040  // Similar to the CoreNEURON solution, we replace all arguments with the
3041  // C++ string that implements them, i.e. `_netcon_data[{}]`. The arguments
3042  // are positional and thus simply numbered through.
3043  const auto& args = node.get_parameters();
3044  RenameVisitor v;
3045  const auto& statement_block = node.get_statement_block();
3046  for (size_t i_arg = 0; i_arg < args.size(); ++i_arg) {
3047  auto old_name = args[i_arg]->get_node_name();
3048  auto new_name = fmt::format("_netcon_data[{}]", i_arg);
3049  v.set(old_name, new_name);
3050  statement_block->accept(v);
3051  }
3052 
3053  auto dparam_it =
3054  std::find_if(info.semantics.begin(), info.semantics.end(), [](const IndexSemantics& a) {
3055  return a.name == naming::FOR_NETCON_SEMANTIC;
3056  });
3057  if (dparam_it == info.semantics.end()) {
3058  throw std::runtime_error("Couldn't find `fornetcon` variable.");
3059  }
3060 
3061  int dparam_index = dparam_it->index;
3062  auto netcon_var = get_name(codegen_int_variables[dparam_index]);
3063 
3064  // This is called from `print_statement_block` which pre-indents the
3065  // current line. Hence `add_text` only.
3066  printer->add_text("double ** _fornetcon_data;");
3067  printer->add_newline();
3068 
3069  printer->fmt_line("int _n_netcons = _nrn_netcon_args({}, &_fornetcon_data);",
3070  get_variable_name(netcon_var, false));
3071 
3072  printer->push_block("for (size_t _i = 0; _i < _n_netcons; ++_i)");
3073  printer->add_line("double * _netcon_data = _fornetcon_data[_i];");
3074  print_statement_block(*statement_block, false, false);
3075  printer->pop_block();
3076 }
3077 
3079  printer->add_line("_NMODLMUTEXLOCK");
3080  printer->add_indent();
3081  node.get_expression()->accept(*this);
3082  printer->add_text(";");
3083  printer->add_line("_NMODLMUTEXUNLOCK");
3084 }
3085 
3086 
3087 } // namespace codegen
3088 } // namespace nmodl
Auto generated AST classes declaration.
Base class for all block scoped nodes.
Definition: block.hpp:41
virtual const ArgumentVector & get_parameters() const
Definition: block.hpp:50
Represents a BREAKPOINT block in NMODL.
Represents a INITIAL block in the NMODL.
Represent LONGITUDINAL_DIFFUSION statement in NMODL.
Definition: lon_diffuse.hpp:39
Extracts information required for LONGITUDINAL_DIFFUSION for each KINETIC block.
const ArgumentVector & get_parameters() const noexcept override
Getter for member variable NetReceiveBlock::parameters.
Base class for all AST node.
Definition: node.hpp:40
Represents block encapsulating list of statements.
Represents a C code block.
Definition: verbatim.hpp:38
Represent WATCH statement in NMODL.
std::vector< std::tuple< std::string, std::string, std::string, std::string > > ParamVector
A vector of parameters represented by a 4-tuple of strings:
std::shared_ptr< symtab::Symbol > SymbolType
std::string table_thread_function_name() const
Name of the threaded table checking function.
void print_nrn_cur() override
Print nrn_cur / current update function definition.
void print_net_receive()
Print net_receive call-back.
void print_mechanism_global_var_structure(bool print_initializers) override
Print the structure that wraps all global variables used in the NMODL.
void print_nrn_destructor() override
Print nrn_destructor function definition.
std::pair< ParamVector, ParamVector > function_table_parameters(const ast::FunctionTableBlock &) override
Parameters of the function itself "{}" and "table_{}".
void print_nrn_current(const ast::BreakpointBlock &node) override
Print the nrn_current kernel.
void print_hoc_py_wrapper_call_impl(const ast::Block *function_or_procedure_block, InterpreterWrapper wrapper_type)
Print the code that calls the impl from the HOC/Py wrapper.
void print_cvode_tolerances()
Print the CVODE function for setup of tolerances.
void print_node_data_structure(bool print_initializers)
Print the structure that wraps all node variables required for the NMODL.
void print_mechanism_range_var_structure(bool print_initializers) override
Print the structure that wraps all range and int variables required for the NMODL.
void print_nrn_state() override
Print nrn_state / state update function definition.
void print_thread_variables_structure(bool print_initializers)
Print the data structure used to access thread variables.
void print_check_table_entrypoint()
Print all check_* function declarations.
std::string global_variable_name(const SymbolType &symbol, bool use_instance=true) const override
Determine the variable name for a global variable given its symbol.
void print_longitudinal_diffusion_callbacks()
Prints the callbacks required for LONGITUDINAL_DIFFUSION.
ParamVector ldifusfunc3_parameters() const
Parameters for what NEURON calls ldifusfunc3_t.
void print_thread_memory_callbacks()
Print thread variable (de-)initialization functions.
std::string get_variable_name(const std::string &name, bool use_instance=true) const override
Determine the C++ string to replace variable names with.
void visit_watch_statement(const ast::WatchStatement &node) override
TODO: Edit for NEURON.
ParamVector cvode_setup_parameters()
Get the parameters for functions that setup (initialize) CVODE.
void print_nrn_cur_kernel(const ast::BreakpointBlock &node) override
Print main body of nrn_cur function.
const ParamVector external_method_parameters(bool table=false) noexcept override
Parameters for functions in generated code that are called back from external code.
std::string register_mechanism_arguments() const override
Arguments for register_mech or point_register_mech function.
void visit_protect_statement(const ast::ProtectStatement &node) override
visit node of type ast::ProtectStatement
void print_nrn_cur_non_conductance_kernel() override
Print the nrn_cur kernel without NMODL conductance keyword provisions.
void print_hoc_py_wrapper(const ast::Block *function_or_procedure_block, InterpreterWrapper wrapper_type)
Print the wrapper for calling FUNCION/PROCEDURES from HOC/Py.
void print_make_node_data() const
Print make_*_node_data.
void visit_for_netcon(const ast::ForNetcon &node) override
visit node of type ast::ForNetcon
void print_global_variables_for_hoc() override
Print byte arrays that register scalar and vector variables for hoc interface.
void print_neuron_includes()
Print includes from NEURON.
void print_entrypoint_setup_code_from_memb_list()
Prints setup code for entrypoints from NEURON.
void print_codegen_routines_nothing()
SUFFIX nothing is special.
ParamVector cvode_update_parameters()
Get the parameters for functions that update state at given timestep in CVODE.
bool optimize_ion_variable_copies() const override
Check if ion variable copies should be avoided.
std::string thread_variable_name(const ThreadVariableInfo &var_info, bool use_instance=true) const
Determine the C++ string to print for thread variables.
void print_global_param_default_values()
Print global struct with default value of RANGE PARAMETERs.
void print_standard_includes() override
Print standard C/C++ includes.
void print_global_function_common_code(BlockType type, const std::string &function_name="") override
Print common code for global functions like nrn_init, nrn_cur and nrn_state.
void print_mechanism_register_regular()
Function body for anything not SUFFIX nothing.
void print_nrn_alloc() override
Print nrn_alloc function definition.
std::string simulator_name() override
Name of the simulator the code was generated for.
void visit_longitudinal_diffusion_block(const ast::LongitudinalDiffusionBlock &node) override
visit node of type ast::LongitudinalDiffusionBlock
void print_cvode_count()
Print the CVODE function returning the number of ODEs to solve.
ParamVector threadargs_parameters()
The parameters for the four macros _threadargs*_.
void print_function_prototypes() override
Print function and procedures prototype declaration.
void print_setdata_functions()
Print NEURON functions related to setting global variables of the mechanism.
std::string backend_name() const override
Name of the code generation backend.
void print_function_or_procedure(const ast::Block &node, const std::string &name, const std::unordered_set< CppObjectSpecifier > &specifiers={ CppObjectSpecifier::Inline}) override
Print nmodl function or procedure (common code)
int position_of_int_var(const std::string &name) const override
Determine the position in the data array for a given int variable.
void print_function_definitions()
Print function and procedures prototype definitions.
void print_entrypoint_setup_code_from_prop()
Prints setup code for entrypoints NEURON.
void print_compute_functions() override
Print all compute functions for every backend.
void print_nrn_constructor() override
Print nrn_constructor function definition.
void print_cvode_definitions()
Print all callbacks for CVODE.
void print_mechanism_register_nothing()
Function body for SUFFIX nothing.
void print_sdlists_init(bool print_initializers) override
void print_verbatim_cleanup(const std::vector< std::string > &macros_defined)
Print #undefs to erase all compatibility macros.
ParamVector internal_method_parameters() override
Parameters for internally defined functions.
std::string get_pointer_name(const std::string &name) const
Determine the C++ string to replace pointer names with.
void print_codegen_routines() override
Print entry point to code generation.
void print_global_macros()
Print NEURON global variable macros.
std::string nrn_thread_internal_arguments() override
Arguments for "_threadargs_" macro in neuron implementation.
std::string py_function_signature(const std::string &function_or_procedure_name) const
Get the signature of the npy <func_or_proc_name> function.
void print_headers_include() override
Print all includes.
void print_macro_definitions()
Print all NEURON macros.
void print_v_unused() const override
Set v_unused (voltage) for NRN_PRCELLSTATE feature.
void visit_verbatim(const ast::Verbatim &node) override
visit node of type ast::Verbatim
ParamVector functor_params() override
The parameters of the Newton solver "functor".
std::string namespace_name() override
Name of "our" namespace.
void print_mechanism_variables_macros()
Print mechanism variables' related macros.
void print_make_instance() const
Print make_*_instance.
void add_variable_tqitem(std::vector< IndexVariableInfo > &variables) override
Add the variable tqitem during get_int_variables.
void print_point_process_function_definitions()
Print POINT_PROCESS related functions Wrap external NEURON functions related to POINT_PROCESS mechani...
void print_fast_imem_calculation() override
Print fast membrane current calculation code.
void print_mechanism_register() override
Print the mechanism registration function.
void print_initial_block(const ast::InitialBlock *node)
Print the initial block.
void print_g_unused() const override
Set g_unused (conductance) for NRN_PRCELLSTATE feature.
std::string py_function_name(const std::string &function_or_procedure_name) const
In non POINT_PROCESS mechanisms all functions and procedures need a py <func_or_proc_name> to be avai...
void print_nrn_jacob()
Print nrn_jacob function definition.
void add_variable_point_process(std::vector< IndexVariableInfo > &variables) override
Add the variable point_process during get_int_variables.
void print_neuron_global_variable_declarations()
Print extern declarations for neuron global variables.
void print_cvode_update(const std::string &name, const ast::StatementBlock &block)
Print the CVODE update function name contained in block.
std::string float_variable_name(const SymbolType &symbol, bool use_instance) const override
Determine the name of a float variable given its symbol.
std::string hoc_function_name(const std::string &function_or_procedure_name) const
All functions and procedures need a hoc <func_or_proc_name> to be available to the HOC interpreter.
void print_nrn_init(bool skip_init_check=true)
Print the nrn_init function definition.
std::string int_variable_name(const IndexVariableInfo &symbol, const std::string &name, bool use_instance) const override
Determine the name of an int variable given its symbol.
void print_net_init()
Print NET_RECEIVE{ INITIAL{ ...
void visit_lon_diffuse(const ast::LonDiffuse &node) override
visit node of type ast::LonDiffuse
std::vector< std::string > print_verbatim_setup(const ast::Verbatim &node, const std::string &verbatim)
Print compatibility macros required for VERBATIM blocks.
void print_net_send_call(const ast::FunctionCall &node) override
Print call to net_send.
const std::string external_method_arguments() noexcept override
Arguments for external functions called from generated code.
void print_net_event_call(const ast::FunctionCall &node) override
Print call to net_event.
ParamVector internalthreadargs_parameters()
The parameters for the four macros _internalthreadargs*_.
void print_global_var_external_access()
Print functions for EXTERNAL use.
std::string hoc_py_wrapper_signature(const ast::Block *function_or_procedure_block, InterpreterWrapper wrapper_type)
Return the wrapper signature.
void print_codegen_routines_regular()
Anything not SUFFIX nothing.
void print_nrn_cur_conductance_kernel(const ast::BreakpointBlock &node) override
Print the nrn_cur kernel with NMODL conductance keyword provisions.
void print_hoc_py_wrapper_setup(const ast::Block *function_or_procedure_block, InterpreterWrapper wrapper_type)
Print the setup code for HOC/Py wrapper.
void append_conc_write_statements(std::vector< ShadowUseStatement > &statements, const Ion &ion, const std::string &concentration) override
Generate Function call statement for nrn_wrote_conc.
void print_data_structures(bool print_initializers) override
Print all classes.
std::string process_verbatim_text(const std::string &verbatim)
ParamVector ldifusfunc1_parameters() const
Parameters for what NEURON calls ldifusfunc1_t.
void print_function_table_call(const ast::FunctionCall &node) override
Print special code when calling FUNCTION_TABLEs.
void print_net_move_call(const ast::FunctionCall &node) override
Print call to net_move.
std::string nrn_thread_arguments() const override
Arguments for "_threadargs_" macro in neuron implementation.
std::string internal_method_arguments() override
Arguments for functions that are defined and used internally.
void print_function_procedure_helper(const ast::Block &node) override
Common helper function to help printing function or procedure blocks.
void print_cvode_setup(const std::string &setup_name, const std::string &update_name)
Print the CVODE setup function setup_name that calls the CVODE update function update_name.
std::string hoc_function_signature(const std::string &function_or_procedure_name) const
Get the signature of the hoc <func_or_proc_name> function.
int position_of_float_var(const std::string &name) const override
Determine the position in the data array for a given float variable.
Class that binds all pieces together for parsing C verbatim blocks.
Definition: c11_driver.hpp:37
Represent symbol table for a NMODL block.
std::vector< std::shared_ptr< Symbol > > get_variables(syminfo::NmodlType with=syminfo::NmodlType::empty, syminfo::NmodlType without=syminfo::NmodlType::empty) const
get variables
bool global_scope() const noexcept
SymbolTable * get_parent_table() const noexcept
Blindly rename given variable to new name
Check if variable is used in given block.
bool variable_used(const ast::Node &node, std::string name)
Visitor for printing C++ code compatible with legacy api of CoreNEURON
Visitor for printing C++ code compatible with legacy api of NEURON
Implement utility functions for codegen visitors.
#define area
Definition: md1redef.h:12
#define v
Definition: md1redef.h:11
#define i
Definition: md1redef.h:19
int thread_data_index
Definition: nocpout.cpp:110
virtual bool is_function_block() const noexcept
Check if the ast node is an instance of ast::FunctionBlock.
Definition: ast.cpp:142
virtual std::string get_node_name() const
Return name of of the node.
Definition: ast.cpp:28
int get_index_from_name(const std::vector< T > &variables, const std::string &name)
BlockType
Helper to represent various block types.
std::string get_name(ast::Ast const *sym)
@ State
derivative block
@ Equation
breakpoint block
static bool ends_with(const std::string &haystack, const std::string &needle)
Check if haystack ends with needle.
static std::vector< std::string > split_string(const std::string &text, char delimiter)
Split a text in a list of words, using a given delimiter character.
std::string join_arguments(const std::string &lhs, const std::string &rhs)
Joint two (list of) arguments.
static bool starts_with(const std::string &haystack, const std::string &needle)
Check if haystack starts with needle.
double(* func)(double)
Definition: hoc_init.cpp:85
double var(InputIterator begin, InputIterator end)
Definition: ivocvect.h:108
#define rhs
Definition: lineq.h:6
const char * name
Definition: init.cpp:16
std::string to_string(EnumT e, const std::array< std::pair< EnumT, std::string_view >, N > &mapping, const std::string_view enum_name)
Converts an enum value to its corresponding string representation.
Definition: nrnreport.hpp:102
static constexpr char POINT_PROCESS_VARIABLE[]
inbuilt neuron variable for point process
static constexpr char NRN_ALLOC_METHOD[]
nrn_alloc method in generated code
static constexpr char AREA_VARIABLE[]
similar to node_area but user can explicitly declare it as area
static constexpr char CVODE_COUNT_NAME[]
name of CVODE method for counting # of ODEs
static constexpr char NRN_CONSTRUCTOR_METHOD[]
nrn_constructor method in generated code
static constexpr char INST_GLOBAL_MEMBER[]
instance struct member pointing to the global variable structure
static constexpr char NRN_INIT_METHOD[]
nrn_init method in generated code
static constexpr char CVODE_SETUP_TOLERANCES_NAME[]
name of CVODE method for setting up tolerances
static constexpr char CONDUCTANCE_UNUSED_VARIABLE[]
range variable when conductance is not used (for vectorized model)
static constexpr char THREAD_ARGS_COMMA[]
verbatim name of the variable for nrn thread arguments, sometimes with trailing comma
static constexpr char THREAD_ARGS[]
verbatim name of the variable for nrn thread arguments
static constexpr char NRN_STATE_METHOD[]
nrn_state method in generated code
static constexpr char USE_TABLE_VARIABLE[]
global variable to indicate if table is used
static constexpr char NTHREAD_DT_VARIABLE[]
dt variable in neuron thread structure
static constexpr char CONDUCTANCE_VARIABLE[]
range variable for conductance
static constexpr char INTERNAL_THREAD_ARGS[]
variation of _threadargs_ for "internal" functions.
static constexpr char NTHREAD_T_VARIABLE[]
t variable in neuron thread structure
static constexpr char FOR_NETCON_SEMANTIC[]
semantic type for for_netcon statement
static constexpr char NRN_CUR_METHOD[]
nrn_cur method in generated code
static constexpr char INTERNAL_THREAD_ARGS_COMMA[]
variation of _threadargs_ for "internal" functions, with comma (maybe).
static constexpr char CVODE_SETUP_NON_STIFF_NAME[]
name of CVODE method for setting up non-stiff systems
static constexpr char CVODE_SETUP_STIFF_NAME[]
name of CVODE method for setting up stiff systems
static constexpr char CVODE_UPDATE_NON_STIFF_NAME[]
name of CVODE method for updating non-stiff systems
static constexpr char NRN_DESTRUCTOR_METHOD[]
nrn_destructor method in generated code
static constexpr char DIAM_VARIABLE[]
inbuilt neuron variable for diameter of the compartment
static constexpr char TQITEM_VARIABLE[]
inbuilt neuron variable for tqitem process
static constexpr char RANDOM_SEMANTIC[]
semantic type for RANDOM variable
static constexpr char ION_VARNAME_PREFIX[]
prefix for ion variable
static constexpr char CVODE_UPDATE_STIFF_NAME[]
name of CVODE method for updating stiff systems
static constexpr char POINTER_SEMANTIC[]
semantic type for pointer variable
static constexpr char NRN_JACOB_METHOD[]
nrn_jacob method in generated code
static constexpr char CVODE_VARIABLE_NAME[]
name of the CVODE variable (can be arbitrary)
static constexpr char VOLTAGE_UNUSED_VARIABLE[]
range variable for voltage when unused (for vectorized model)
static constexpr char NODE_AREA_VARIABLE[]
inbuilt neuron variable for area of the compartment
static void rename_net_receive_arguments(const ast::NetReceiveBlock &net_receive_node, const ast::Node &node)
Rename arguments to NET_RECEIVE block with corresponding pointer variable.
std::unordered_map< std::string, std::string > get_nonglobal_local_variable_names(const symtab::SymbolTable &symtab)
Map of the non-(global/top-local) LOCAL variables.
InterpreterWrapper
Enum to switch between HOC and Python wrappers for functions and procedures defined in mechanisms.
NmodlType
NMODL variable properties.
encapsulates code generation backend implementations
Definition: ast_common.hpp:26
static List * info
Symbol * breakpoint_current(Symbol *s)
Definition: nocpout.cpp:3288
void function_table(Symbol *s, Item *qpar1, Item *qpar2, Item *qb1, Item *qb2)
Definition: parsact.cpp:940
static Node * node(Object *)
Definition: netcvode.cpp:291
return status
size_t p
short type
Definition: cabvars.h:10
static List * current
Definition: nrnunit.cpp:13
static int point_process
Definition: nrnunit.cpp:12
#define text
Definition: plot.cpp:60
Auto generated AST classes declaration.
Blindly rename given variable to new name
static uint32_t value
Definition: scoprand.cpp:25
Implement string manipulation functions.
Represent semantic information for index variable.
Helper to represent information about index/int variables.
bool is_integer
if this is an integer (e.g.
bool is_index
if this is pure index (e.g.
Represent ions used in mod file.
std::string rev_potential_pointer_name() const
std::string name
name of the ion
std::string extra_conc_pointer_name() const
std::string intra_conc_pointer_name() const
Represents ion write statement during code generation.
const std::shared_ptr< symtab::Symbol > symbol
size_t offset
The global variables ahead of this one require offset doubles to store.
Definition: units.cpp:71
nmodl::parser::UnitDriver driver
Definition: parser.cpp:28
Check if variable is used in given block.
Utility functions for visitors implementation.