File: programmatic.rst

package info (click to toggle)
neuron 8.2.6-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 34,760 kB
  • sloc: cpp: 149,571; python: 58,465; ansic: 50,329; sh: 3,510; xml: 213; pascal: 51; makefile: 35; sed: 5
file content (493 lines) | stat: -rw-r--r-- 14,174 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
Programmatic Simulation Control
===============================

See also:

.. toctree:: :maxdepth: 1

    cvode.rst
    batch.rst
    savstate.rst
    bbsavestate.rst
    sessionsave.rst

Functions
---------

.. hoc:function:: initnrn


    Syntax:
        ``initnrn()``


    Description:
        Initialize ``t, dt, clamp_resist``, and ``celsius`` to the values 
        they had when the program was first run. 
         
        Note that in this 
        version ``Ra`` is no longer a global variable but a section variable 
        like *L* and *rallbranch*. Thus ``Ra`` can be different for different 
        sections.  In order to set ``Ra`` to a constant value, use: 
         
        ``forall Ra=...`` 

    .. warning::
        Not very useful. No way to completely restart neuron except to :hoc:func:`quit` and
        re-load. 


----



.. hoc:function:: fadvance


    Syntax:
        ``fadvance()``


    Description:
        Integrate all section equations over the interval :hoc:data:`dt` .
        The value of :hoc:data:`t` is incremented by dt.
        The default method is first order implicit but may be changed to 
        Crank-Nicolson by changing :hoc:data:`secondorder` = 2.
         
        fadvance integrates the equation over the dt step by 
        calling all the BREAKPOINT blocks of models at t+dt/2 twice with 
        v+.001 and v in order to compute the current and conductance to form 
        the matrix conductance*voltage = current. 
        This matrix is then solved for v(t+dt). 
        (if secondorder == 2 the ionic currents are adjusted to be second order 
        correct. If secondorder == 1 the ionic currents are not adjusted but 
        the voltages are second order correct) 
        Lastly the SOLVE statement within the BREAKPOINT block of models is 
        executed with t+dt and the new values of v in order to integrate those 
        states (from new t-.5dt to new t+.5dt). 

         

----



.. hoc:function:: finitialize


    Syntax:
        ``finitialize()``

        ``finitialize(v)``


    Description:
        Call the INITIAL block for all mechanisms and point processes 
        inserted in the sections. 
        If the optional argument is present then all voltages of all sections 
        are initialized to *v*. 
        :hoc:data:`t` is set to 0.
         
        The order of principal actions during an finitialize call is:
        
        -   Type 3 FInitializeHandler statements executed. 
        -   Make sure internal structures needed by integration methods are consistent 
             with the current biophysical spec. 
        -   t = 0 
        -   Clear the event queue. 
        -   Random.play values assigned to variables. 
        -   Vector.play at t=0 values assigned to variables. 
        -   All v = arg if the arg is present. 
        -   Type 0 FInitializeHandler statements executed. 
        -   All mechanism BEFORE INITIAL blocks are called. 
        -   All mechanism INITIAL blocks called. 
               Mechanisms that WRITE concentrations are after ion mechanisms and 
               before mechanisms that READ concentrations. 
        -   LinearMechanism states are initialized 
        -   INITIAL blocks inside NETRECEIVE blocks are called. 
        -   All mechanism AFTER INITIAL blocks are called. 
        -   Type 1 FInitializeHandler statements executed. 
        -   The INITIAL block net_send(0, flag) events are delivered. 
        -   Effectively a call to CVode.re_init or fcurrent(), whichever appropriate. 
        -   Various record functions at t=0. e.g. CVode.record, Vector.record  
        -   Type 2 FInitializeHandler statements executed. 
             


    .. seealso::
        :hoc:class:`FInitializeHandler`, :ref:`hoc_runcontrol_Init`, :hoc:meth:`CVode.re_init`, :hoc:func:`fcurrent`, :hoc:func:`frecord_init`

         

----



.. hoc:function:: frecord_init


    Syntax:
        ``frecord_init()``


    Description:
        Initializes the Vectors which are recording variables. i.e. resize to 0 and 
        append the current values of the variables.  This is done at the end 
        of an :hoc:func:`finitialize` call but needs to be done again to complete initialization
        if the user changes states or assigned variables that are being recorded.. 

    .. seealso::
        :hoc:meth:`Vector.record`, :ref:`hoc_runcontrol_Init`

----





.. hoc:function:: fcurrent


    Syntax:
        ``fcurrent()``


    Description:
        Make all assigned variables (currents, conductances, etc) 
        consistent with the values of the states. Useful in combination 
        with :hoc:func:`finitialize`.

    Example:

        .. code-block::
            none

            create soma 
            access soma 
            insert hh 
            print "default el_hh = ", el_hh 
            // set el_hh so that the steady state is exactly -70 mV 
            finitialize(-70) // sets v to -70 and m,h,n to corresponding steady state values 
             
            fcurrent()	// set all assigned variables consistent with states 
             
            // use current balance: 0 = ina + ik + gl_hh*(v - el_hh)		 
            el_hh = (ina + ik + gl_hh*v)/gl_hh 
             
            print "-70 mV steady state el_hh = ", el_hh 
            fcurrent()	// recalculate currents (il_hh) 


         

----



.. hoc:function:: fmatrix


    Syntax:
        ``fmatrix()``

        ``section {value = fmatrix(x, index)}``


    Description:
        No args: print the jacobian matrix for the tree structure in a particularly 
        confusing way. for debugging only. 
         
        With args, return the matrix element associated with the integer index 
        in the row corresponding to the currently accessed 
        section at position x. The index 1...4 is associated with: 
        The coeeficient for the effect of this locations voltage on current balance at the parent location, 
        The coeeficient for the effect of this locations voltage on current balance at this location, 
        The coeeficient for the effect of the parent locations voltage on current balance at this location, 
        The right hand side of the matrix equation for this location. These are the 
        values of NODEA, NODED NODEB, and NODERHS respectively in 
        nrn/src/nrnoc/section.h . The matrix elements are properly setup on return 
        from a call to the :hoc:func:`fcurrent` function. For the fixed step method
        :hoc:func:`fadvance` modifies NODED and NODERHS
        but leaves NODEA and NODEB unchanged. 

         
----

.. hoc:data:: secondorder


    Syntax:
        ``secondorder``


    Description:
        This is a global variable which specifies the time integration method. 


        =0 
            default fully implicit backward euler. Very numerically stable. 
            gives steady state in one step when *dt=1e10*. Numerical errors 
            are proportional to :hoc:data:`dt`.

        =1 
            crank-nicholson Can give large (but damped) numerical error 
            oscillations. For small :hoc:data:`dt` the numerical errors are proportional
            to ``dt^2``. Cannot be used with voltage clamps. Ionic currents 
            are first order correct. Channel conductances are second order 
            correct when plotted at ``t+dt/2`` 

        =2 
            crank-nicholson like 1 but in addition Ion currents (*ina*, *ik*, 
            etc) are fixed up so that they are second order correct when 
            plotted at ``t-dt/2`` 


         

----



.. hoc:data:: t


    Syntax:
        ``t``


    Description:
        The global time variable. 

         

----



.. hoc:data:: dt


    Syntax:
        ``dt``


    Description:
        The integration interval for :hoc:func:`fadvance`.
         
        When using the default implicit integration method (:hoc:data:`secondorder` = 0)
        there is no upper limit on dt for numerical stability and in fact for 
        passive models it is often convenient to use dt=1.9 to obtain the 
        steady state in a single time step. 
         
        dt can be changed by the user at any time during a simulation. However, 
        some inserted mechanisms may use tables which depend on the value of dt 
        which will be automatically recomputed. In this situation, the tables 
        are not useful and should be bypassed by setting the appropriate 
        usetable_suffix global variables to 0. 

         

----



.. hoc:data:: clamp_resist


    Syntax:
        ``clamp_resist``


    Description:
        Obsolete, used by fclamp. 

         

----



.. hoc:data:: celsius


    Syntax:
        ``celsius = 6.3``


    Description:
        Temperature in degrees centigrade. 
         
        Generally, rate function tables (eg. used by the hh mechanism) 
        depend on temperature and will automatically be re-computed 
        whenever celsius changes. 

         

----



.. hoc:data:: stoprun


    Syntax:
        ``stoprun``


    Description:
        A flag which is watched by :hoc:func:`fit_praxis`, :hoc:class:`CVode`, and other procedures
        during a run or family of runs. 
        When stoprun==1 they will immediately return without completing 
        normally. This allows safe stopping in the middle of a long run. Procedures 
        that do multiple runs should check stoprun after each run and exit 
        gracefully. The :hoc:meth:`RunControl.Stop` of the RunControl GUI sets this variable.
        It is cleared at the beginning of a run or when continuing a run. 


----

.. hoc:function:: checkpoint

    Syntax:
        :samp:`checkpoint("{filename}")`

    Description:
        saves the current state of the system in a portable file to 
        allow one to take up where you left off -- possibly on another 
        machine. Returning to this state is accomplished by running the 
        program with the checkpoint file as the first argument. 
        If the checkpoint file is inconsistent with the executable the 
        program prints an error message and exits. 
         
        At this time many portions of the computer state are left out of the 
        checkpoint file, i.e. it is not as complete as a core dump. 
        Some things that ARE included are: 
        all interpreter symbols with definitions and values, 
        all hoc instructions, 
        all neuron state/parameters with mechanisms. 
        Many aspects of the GUI are not included. 
         
    .. warning::
        There is not enough implementation at this time to make this 
        facility useful. Use the :hoc:class:`SaveState` class instead.



         
----


.. _hoc_finithnd:

         
FInitializeHandler
------------------



.. hoc:class:: FInitializeHandler


    Syntax:
        ``fih = new FInitializeHandler("stmt", [obj])``

        ``fih = new FInitializeHandler(type, "stmt", [obj])``


    Description:
        Install an initialization handler statement to be called during a call to 
        :hoc:func:`finitialize`. The default type is 1. The
        statement will be executed at the top level of the interpreter 
        or else in the context of the optional obj arg. 
         
        Type 0 handlers are called before the mechanism INITIAL blocks. 
         
        Type 1 handlers are called after the mechanism INITIAL blocks. 
        This is the best place to change state values. 
         
        Type 2 handlers are called just before return from finitialize. 
        This is the best place to record values at t=0. 
         
        Type 3 handlers are called at the beginning of finitialize. 
        At this point it is allowed to change the structure of the model. 
         
        See :hoc:func:`finitialize` for more details about the order of initialization processes
        within that function. 
         
        This class helps alleviate the administrative problems of maintaining variations 
        of the proc :ref:`hoc_RunControl_Init`.

    Example:

        .. code-block::
            none

            // specify an example model 
            load_file("nrngui.hoc") 
            create a, b 
            access a 
            forall insert hh 
             
            objref fih[3] 
            fih[0] = new FInitializeHandler(0, "fi0()") 
            fih[1] = new FInitializeHandler(1, "fi1()") 
            fih[2] = new FInitializeHandler(2, "fi2()") 
             
            proc fi0() { 
            	print "fi0() called after v set but before INITIAL blocks" 
            	printf("  a.v=%g a.m_hh=%g\n", a.v, a.m_hh) 
            	a.v = 10 
            } 
             
            proc fi1() { 
            	print "fi1() called after INITIAL blocks but before BREAKPOINT blocks" 
            	print "     or variable step initialization." 
            	print "     Good place to change any states." 
            	printf("  a.v=%g a.m_hh=%g\n", a.v, a.m_hh) 
            	printf("  b.v=%g b.m_hh=%g\n", b.v, b.m_hh) 
            	b.v = 10 
            } 
             
            proc fi2() { 
            	print "fi2() called after everything initialized. Just before return" 
            	print "     from finitialize." 
            	print "     Good place to record or plot initial values" 
            	printf("  a.v=%g a.m_hh=%g\n", a.v, a.m_hh) 
            	printf("  b.v=%g b.m_hh=%g\n", b.v, b.m_hh) 
            } 
             
            begintemplate Test 
            objref fih, this 
            proc init() { 
            	fih = new FInitializeHandler("p()", this) 
            } 
            proc p() { 
            	printf("inside %s.p()\n", this) 
            } 
            endtemplate Test 
             
            objref test 
            test = new Test() 
             
            stdinit() 
            fih[0].allprint() 


         

----



.. hoc:method:: FInitializeHandler.allprint


    Syntax:
        ``fih.allprint()``


    Description:
        Prints all the FInitializeHandler statements along with their object context 
        in the order they will be executed during an :hoc:func:`finitialize` call.