OB.DAAC Logo
NASA Logo
Ocean Color Science Software

ocssw V2022
ffnet.f
Go to the documentation of this file.
1 c#######################################################################
2 c# Copyright (C) 2006 by Marek Wojciechowski
3 c# <mwojc@p.lodz.pl>
4 c#
5 c# Distributed under the terms of the GNU General Public License (GPL)
6 c# http://www.gnu.org/copyleft/gpl.html
7 c#######################################################################
8 c
9 cc
10 ccc
11 cccc
12 ccccc BASIC FFNN PROPAGATION ROUTINES
13 cccc
14 ccc
15 cc
16 c
17 c************************************************************************
18  subroutine prop(x, conec, n, units, u)
19 c************************************************************************
20 c
21 c.....Gets conec and units with input already set
22 c.....and calculates all activations.
23 c.....Identity input and sigmoid activation function for other units
24 c.....is assumed
25 c
26  implicit none
27 c.....variables
28  integer n, u, conec(n,2)
29  double precision x(n), units(u)
30 c.....helper variables
31  integer src, trg, ctrg, xn
32 c.....f2py statements
33 cf2py intent(in, out) units
34 
35 c.....propagate signals with sigmoid activation function
36  if (n.gt.0) then
37  ctrg = conec(1,2)
38  units(ctrg) = 0.
39  do xn=1,n
40  src = conec(xn,1)
41  trg = conec(xn,2)
42  !if next unit
43  if (trg.ne.ctrg) then
44  units(ctrg) = 1./(1.+exp(-units(ctrg)))
45  ctrg = trg
46  if (src.eq.0) then !handle bias
47  units(ctrg) = x(xn)
48  else
49  units(ctrg) = units(src) * x(xn)
50  endif
51  else
52  if (src.eq.0) then !handle bias
53  units(ctrg) = units(ctrg) + x(xn)
54  else
55  units(ctrg) = units(ctrg) + units(src) * x(xn)
56  endif
57  endif
58  enddo
59  units(ctrg) = 1./(1.+exp(-units(ctrg))) !for last unit
60  endif
61 
62  return
63  end
64 
65 c************************************************************************
66  subroutine sqerror(x, conec, n, units, u, inno, i, outno, o,
67  & Input, Targ, p, sqerr)
68 c************************************************************************
69 c
70 c.....Takes Input and Target patterns and returns sum of squared errors
71 c
72  implicit none
73 c.....variables
74  integer n, u, i, o, p, conec(n,2), inno(i), outno(o)
75  double precision x(n), units(u), Input(p,i), Targ(p,o), sqerr
76 c.....helper variables
77  integer k, pat
78 c.....f2py statements
79 cf2py intent(out) sqerr
80 
81  sqerr = 0.
82 c.....loop over given patterns
83  DO pat = 1,p
84 c.........set input units
85  do k=1,i
86  units(inno(k)) = input(pat,k)
87  enddo
88 c.........propagate signals
89  call prop(x, conec, n, units, u)
90 c.........sum squared errors
91  do k=1,o
92  sqerr = sqerr + (units(outno(k)) - targ(pat,k))**2
93  enddo
94  ENDDO
95  sqerr = 0.5d0*sqerr
96 
97  return
98  end
99 
100 c************************************************************************
101  subroutine grad(x, conec, n, bconecno, bn, units, u, inno, i,
102  & outno, o, Input, Targ, p, xprime)
103 c************************************************************************
104 c
105 c.....Takes conec, bconecno, Input and Target patterns and returns
106 c.....gradient calculated with error backpropagation
107 c
108  implicit none
109 c.....variables
110  integer n, bn, u, i, o, p
111  integer conec(n,2), bconecno(bn), inno(i), outno(o)
112  double precision x(n), units(u), Input(p,i), Targ(p,o)
113  double precision xprime(n), diff(o), bunits(u)
114 c.....helper variables
115  integer k, pat, src, trg, ctrg
116  double precision deriv, cx
117 c.....f2py statements
118 cf2py intent(out) xprime
119 
120 c.....initialize xprime
121  do k=1,n
122  xprime(k) = 0.
123  enddo
124 c.....loop over given patterns
125  DO pat = 1,p
126 c.........propagate input signals
127  do k=1,i
128  units(inno(k)) = input(pat,k)
129  enddo
130  call prop(x, conec, n, units, u)
131 c.........set diffs at network output as back network inputs
132  do k=1,o
133  diff(k) = units(outno(k)) - targ(pat,k)
134  deriv = units(outno(k)) * (1. - units(outno(k))) !ugly
135  bunits(outno(k)) = diff(k) * deriv
136  enddo
137 c.........backpropagate errors
138  if (bn.gt.0) then
139  ctrg = conec(bconecno(1),1)
140  bunits(ctrg) = 0.
141  do k=1,bn
142  src = conec(bconecno(k),2)
143  trg = conec(bconecno(k),1)
144  cx = x(bconecno(k))
145  if (trg.ne.ctrg) then !if next unit
146  deriv = units(ctrg) * (1. - units(ctrg)) !ugly
147  bunits(ctrg) = bunits(ctrg) * deriv
148  ctrg = trg
149  bunits(ctrg) = bunits(src) * cx
150  else
151  bunits(ctrg) = bunits(ctrg) + bunits(src) * cx
152  endif
153  enddo
154  deriv = units(ctrg) * (1 - units(ctrg))
155  bunits(ctrg) = bunits(ctrg) * deriv !for last unit
156  endif
157 c.........add gradient elements to overall xprime
158  do k=1,n
159  src = conec(k,1)
160  trg = conec(k,2)
161  if (src.eq.0) then !handle bias
162  xprime(k) = xprime(k) + bunits(trg)
163  else
164  xprime(k) = xprime(k) + units(src) * bunits(trg)
165  endif
166  enddo
167  ENDDO
168 
169  return
170  end
171 
172 c************************************************************************
173  subroutine recall(x, conec, n, units, u, inno, i, outno, o,
174  & input, output)
175 c************************************************************************
176 c
177 c.....Takes single input pattern and returns network output
178 c
179  implicit none
180 c.....variables
181  integer n, u, i, o, conec(n,2), inno(i), outno(o)
182  double precision x(n), units(u), input(i), output(o)
183 c.....helper variables
184  integer k
185 c.....f2py statements
186 cf2py intent(out) output
187 
188 c.....set input units
189  do k=1,i
190  units(inno(k)) = input(k)
191  enddo
192 c.....propagate signals
193  call prop(x, conec, n, units, u)
194 c.....get output
195  do k=1,o
196  output(k) = units(outno(k))
197  enddo
198 
199  return
200  end
201 
202 c************************************************************************
203  subroutine diff(x, conec, n, dconecno, dn, dconecmk, units, u,
204  & inno, i, outno, o, input, deriv)
205 c************************************************************************
206 c
207 c.....Takes single input pattern and returns network partial derivatives
208 c.....in the form d(output,o)/d(input,i). 'units' contain now activation
209 c.....derivatives
210 c
211  implicit none
212 c.....variables
213  integer n, dn, u, i, o, conec(n,2), dconecno(dn), dconecmk(i+1)
214  integer inno(i), outno(o)
215  double precision x(n), units(u), dunits(u), input(i), deriv(o,i)
216 c.....helper variables
217  integer k, di, trg, src, ctrg, xn
218  double precision dx
219 c.....f2py statements
220 cf2py intent(out) deriv
221 
222 c.....first set inputs for usual and derivative network units
223  do k=1,i
224  units(inno(k)) = input(k)
225  dunits(inno(k)) = 0d0
226  enddo
227 c.....calculate derivatives of activation functions --> units became
228 c.....units derivatives (ugly, usable only for sigmoid function
229 c.....and identity input)
230  call prop(x, conec, n, units, u)
231  do k=1,u
232  units(k) = units(k) * (1d0 - units(k))
233  enddo
234  do k=1,i
235  units(inno(k)) = 1d0
236  enddo
237 c.....prepare output units for derivative network
238  do k=1,o
239  dunits(outno(k)) = 0d0
240  enddo
241 c.....loop over inputs
242  do di=1,i
243 c.........set current input unit derivative as network input
244  dunits(inno(di)) = units(inno(di))
245 c.........propagate signals through derivative network (dunits became
246 c.........net units and units derivatives are now scaling factors)
247  if (dn.gt.0) then
248  ctrg = conec(dconecno(dconecmk(di)+1),2)
249  dunits(ctrg) = 0.
250  do xn=dconecmk(di)+1,dconecmk(di+1)
251  src = conec(dconecno(xn),1)
252  trg = conec(dconecno(xn),2)
253  dx = x(dconecno(xn))
254  if (trg.ne.ctrg) then
255  dunits(ctrg) = dunits(ctrg) * units(ctrg)
256  ctrg = trg
257  dunits(ctrg) = dunits(src) * dx
258  else
259  dunits(ctrg) = dunits(ctrg) + dunits(src) * dx
260  endif
261  enddo
262  dunits(ctrg) = dunits(ctrg) * units(ctrg) !for last unit
263  endif
264 c.........save network outputs (do/di)
265  do k=1,o
266  deriv(k, di) = dunits(outno(k))
267  dunits(outno(k)) = 0d0
268  enddo
269 c.........restore current input
270  dunits(inno(di)) = 0d0
271  enddo
272 
273  RETURN
274  end
275 c
276 cc
277 ccc
278 cccc
279 ccccc EXTENSIONS OF BASIC ROUTINES
280 cccc
281 ccc
282 cc
283 c
284 c************************************************************************
285  subroutine func(x, conec, n, bconecno, bn, units, u, inno, i,
286  & outno, o, Input, Targ, p, sqerr)
287 c************************************************************************
288 c
289 c.....Just calls sqerror, but now the agruments list
290 c.....is compatibile with grad. This compatibility is needed by scipy
291 c.....optimizers.
292 c
293  implicit none
294 c.....variables
295  integer n, bn, u, i, o, p
296  integer conec(n,2), bconecno(bn), inno(i), outno(o)
297  double precision x(n), units(u), Input(p,i), Targ(p,o), sqerr
298 c.....f2py statements
299 cf2py intent(out) sqerr
300 
301  call sqerror(x, conec, n, units, u, inno, i, outno, o,
302  & input, targ, p, sqerr)
303 
304  return
305  end
306 
307 c************************************************************************
308  subroutine pikaiaff(x, ffn, conec, n, units, u, inno, i, outno, o,
309  & Input, Targ, p, bound1, bound2, isqerr)
310 c************************************************************************
311 c
312 c.....Routine for use with pikaia - genetic algorithm based optimizer.
313 c.....Takes Input and Target patterns and returns inverse of
314 c.....sum of quared errors. Note: (bound1, bound2)
315 c.....is constraint range for x.
316 c
317  implicit none
318 c.....variables
319  integer n, ffn, u, i, o, p, conec(n,2), inno(i), outno(o)
320  double precision x(n), x2(n), units(u), Input(p,i), Targ(p,o)
321  double precision bound1, bound2, isqerr
322 c.....f2py statements
323 cf2py intent(out) isqerr
324 
325 c.....first map x vector values from 0,1 to bound1,bound2
326  call vmapa(x, n, 0d0, 1d0, bound1, bound2, x2)
327 c.....now propagate patterns and obtain error
328  call sqerror(x2, conec, n, units, u, inno, i, outno, o,
329  & input, targ, p, isqerr)
330 c.....inverse error
331  isqerr = 1. / isqerr
332 
333  RETURN
334  end
335 
336 c************************************************************************
337  subroutine normcall(x, conec, n, units, u, inno, i, outno, o,
338  & eni, deo, input, output)
339 c************************************************************************
340 c
341 c.....Takes single input pattern and returns network output
342 c.....This have the same functionality as recall but now input and
343 c.....output are normalized inside the function.
344 c.....eni = [ai, bi], eno = [ao, bo] - parameters of linear mapping
345 c
346  implicit none
347 c.....variables
348  integer n, u, i, o, conec(n,2), inno(i), outno(o)
349  double precision x(n), units(u), input(i), output(o)
350  double precision eni(i,2), deo(o,2)
351 c.....f2py statements
352 cf2py intent(out) output, istat
353 
354 c.....set input units
355  call setin(input, inno, i, eni, units, u)
356 c.....propagate signals
357  call prop(x, conec, n, units, u)
358 c.....get output
359  call getout(units, u, outno, o, deo, output)
360 
361  return
362  end
363 c
364 c************************************************************************
365  subroutine normdiff(x, conec, n, dconecno, dn, dconecmk, units,
366  & u, inno, i, outno, o, eni, ded, input, deriv)
367 c************************************************************************
368 c
369 c.....Takes single input pattern and returns network partial derivatives
370 c.....in the form d(output,o)/d(input,i). 'units' contain now activation
371 c.....derivatives
372 c.....This have the same functionality as diff but now input and
373 c.....output are normalized inside function
374 c
375 c.....Solution not very smart, whole diff routine is rewritten here...
376 c
377  implicit none
378 c.....variables
379  integer n, dn, u, i, o, conec(n,2), dconecno(dn), dconecmk(i+1)
380  integer inno(i), outno(o)
381  double precision x(n), units(u), dunits(u), input(i), deriv(o,i)
382  double precision eni(i,2), ded(o,i)
383 c.....helper variables
384  integer k, di, trg, src, ctrg, xn
385  double precision dx
386 c.....f2py statements
387 cf2py intent(out) deriv
388 
389 c.....first set inputs for usual and derivative network units
390  call setin(input, inno, i, eni, units, u)
391 c.....propagate through network
392  call prop(x, conec, n, units, u)
393 c.....calculate derivatives of activation functions --> units became
394 c.....units derivatives (ugly, usable only for sigmoid function
395 c.....and identity input)
396  do k=1,u
397  units(k) = units(k) * (1d0 - units(k))
398  !units(k) = units(k) * (1d0 - units(k)) * (1d0-2d0*units(k))
399  enddo
400  do k=1,i
401  units(inno(k)) = 1d0
402  enddo
403 c.....prepare output units and scaling factors for derivative network
404  do k=1,o
405  dunits(outno(k)) = 0d0
406  enddo
407 c.....loop over inputs
408  do di=1,i
409 c.........set current input unit derivative as network input
410  dunits(inno(di)) = units(inno(di))
411 c.........propagate signals through derivative network (dunits became
412 c.........net units and units derivatives are now scaling factors)
413  if (dn.gt.0) then
414  ctrg = conec(dconecno(dconecmk(di)+1),2)
415  dunits(ctrg) = 0.
416  do xn=dconecmk(di)+1,dconecmk(di+1)
417  src = conec(dconecno(xn),1)
418  trg = conec(dconecno(xn),2)
419  dx = x(dconecno(xn))
420  if (trg.ne.ctrg) then
421  dunits(ctrg) = dunits(ctrg) * units(ctrg)
422  ctrg = trg
423  dunits(ctrg) = dunits(src) * dx
424  else
425  dunits(ctrg) = dunits(ctrg) + dunits(src) * dx
426  endif
427  enddo
428  dunits(ctrg) = dunits(ctrg) * units(ctrg) !for last unit
429  endif
430 c.........save network outputs (do/di)
431  do k=1,o
432  deriv(k, di) = dunits(outno(k))*ded(k,di)
433  dunits(outno(k)) = 0d0
434  enddo
435 c.........restore current input
436  dunits(inno(di)) = 0d0
437  enddo
438 
439  RETURN
440  end
441 c
442 c************************************************************************
443  subroutine normcall2(x, conec, n, units, u, inno, i, outno, o,
444  & eni, deo, input, p, output)
445 c************************************************************************
446 c
447 c.....Calls normcall for an array if inputs and return array of outputs
448 c
449  implicit none
450 c.....variables
451  integer n, u, i, o, conec(n,2), inno(i), outno(o), p
452  double precision x(n), units(u), input(p,i), output(p,o)
453  double precision eni(i,2), deo(o,2)
454 c.....helper variables
455  double precision tmpinp(i), tmpout(o)
456  integer j, k
457 c.....f2py statements
458 cf2py intent(out) output, istat
459 
460 c.....iterate over input set
461  do j=1,p
462  do k=1,i
463  tmpinp(k) = input(j,k)
464  enddo
465  call normcall(x, conec, n, units, u, inno, i, outno, o,
466  & eni, deo, tmpinp, tmpout)
467  do k=1,o
468  output(j,k) = tmpout(k)
469  enddo
470  enddo
471 
472  return
473  end
474 c
475 c************************************************************************
476  subroutine normdiff2(x, conec, n, dconecno, dn, dconecmk, units,u,
477  & inno, i, outno, o, eni, ded, input, p, deriv)
478 c************************************************************************
479 c
480 c.....Calls normdiff for an array if inputs and return array of derivs
481 c
482  implicit none
483 c.....variables
484  integer n, dn, u, i, o, p, conec(n,2), dconecno(dn), dconecmk(i+1)
485  integer inno(i), outno(o)
486  double precision x(n), units(u), input(p,i), deriv(p,o,i)
487  double precision eni(i,2), ded(o,i)
488 c.....helper variables
489  integer j, k, l
490  double precision tmpinp(i), tmpder(o,i)
491 c.....f2py statements
492 cf2py intent(out) deriv
493 
494 c.....iterate over input set
495  do j=1,p
496  do k=1,i
497  tmpinp(k) = input(j,k)
498  enddo
499  call normdiff(x, conec, n, dconecno, dn, dconecmk, units,
500  & u, inno, i, outno, o, eni, ded, tmpinp, tmpder)
501  do k=1,o
502  do l=1,i
503  deriv(j,k,l) = tmpder(k,l)
504  enddo
505  enddo
506  enddo
507 
508  return
509  end
510 c
511 cc
512 ccc
513 cccc
514 ccccc BASIC TRAINING ALGORITHMS
515 cccc
516 ccc
517 cc
518 c
519 c************************************************************************
520  subroutine momentum(x, conec, n, bconecno, bn, units, u, inno, i,
521  & outno, o, Input, Targ, p, eta, moment, maxiter)
522 c************************************************************************
523 c
524 c.....Standard backpropagation training with momentum
525 c
526  implicit none
527 c.....variables
528  integer n, bn, u, i, o, p, maxiter
529  integer conec(n,2), bconecno(bn), inno(i), outno(o)
530  double precision x(n), units(u), Input(p,i), Targ(p,o)
531  double precision xprime(n), update, update0(n), eta, moment
532 c.....helper variables
533  integer j, k
534 c.....f2py statements
535 cf2py intent(in, out) x
536 
537 c.....initialize variables
538  do j=1,n
539  update0(j) = 0d0
540  enddo
541  k=0
542 c.....update maxiter times
543  do while (k.lt.maxiter)
544  call grad(x, conec, n, bconecno, bn, units, u, inno, i,
545  & outno, o, input, targ, p, xprime)
546  do j=1,n
547  update = -eta*xprime(j)
548  x(j) = x(j) + update + moment*update0(j)
549  update0(j) = update
550  enddo
551  k=k+1
552  enddo
553 
554  return
555  end
556 
557 c************************************************************************
558  subroutine rprop(x, conec, n, bconecno, bn, units, u, inno, i,
559  & outno, o, Input, Targ, p,
560  & a, b, mimin, mimax, xmi, maxiter)
561 c************************************************************************
562 c
563 c.....Rprop training algorithm
564 c
565  implicit none
566 c.....variables
567  integer n, bn, u, i, o, p, maxiter
568  integer conec(n,2), bconecno(bn), inno(i), outno(o)
569  double precision x(n), units(u), Input(p,i), Targ(p,o)
570  double precision xprime(n), xprime0(n), xmi(n)
571  double precision a, b, mimax, mimin
572 c.....helper variables
573  integer j, k
574 c.....f2py statements
575 cf2py intent(in, out) x, xmi
576 
577 c.....initialize variables
578  do j=1,n
579  xprime0(j) = 0d0
580  enddo
581  k=0
582 c.....update maxiter times
583  do while (k.lt.maxiter)
584  call grad(x, conec, n, bconecno, bn, units, u, inno, i,
585  & outno, o, input, targ, p, xprime)
586  do j=1,n
587 c.............find mi coefficient
588  if ( xprime(j) * xprime0(j) .gt. 0 ) then
589  xmi(j) = min( a * xmi(j), mimax )
590  elseif ( xprime(j) * xprime0(j) .lt. 0 ) then
591  xmi(j) = max( b * xmi(j), mimin )
592  else
593  xmi(j) = xmi(j)
594  endif
595 c.............update weights and record gradient components
596  x(j) = x(j) - sign( xmi(j), xprime(j) )
597  xprime0(j) = xprime(j)
598  enddo
599  k=k+1
600  enddo
601 
602  return
603  end
604 c
605 cc
606 ccc
607 cccc
608 ccccc HELPER FUNCTIONS AND ROUTINES
609 cccc
610 ccc
611 cc
612 c
613 c************************************************************************
614  subroutine setin(input, inno, i, eni, units, u)
615 c************************************************************************
616 c
617 c.....normalize and set input units
618 c
619  implicit none
620 c.....variables
621  integer i, u, inno(i)
622  double precision input(i), units(u), eni(i,2)
623 c.....helper variables
624  integer k
625 c.....f2py statements
626 cf2py intent(in,out) units
627 
628  do k=1,i
629  units(inno(k)) = eni(k,1) * input(k) + eni(k,2)
630  enddo
631 
632  return
633  end
634 
635 c************************************************************************
636  subroutine getout(units, u, outno, o, deo, output)
637 c************************************************************************
638 c
639 c.....get and denormalize output units
640 c
641  implicit none
642 c.....variables
643  integer o, u, outno(o)
644  double precision output(o), units(u), deo(o,2)
645 c.....helper variables
646  integer k
647 c.....f2py statements
648 cf2py intent(out) output
649 
650  do k=1,o
651  output(k) = deo(k,1) * units(outno(k)) + deo(k,2)
652  enddo
653 
654  return
655  end
656 
657 c************************************************************************
658  function mapa(f, a, b, c, d)
659 c************************************************************************
660 c
661 c.....Linear map of f from range (a,b) to range (c,d)
662 c
663  implicit none
664 c.....variables
665  double precision a, b, c, d, f, mapa
666 c.....helper variables
667  double precision t
668 
669 c.....map vector (no check of bounds...)
670  t = ( d - c ) / ( b - a )
671  mapa = c + ( f - a ) * t
672 
673  RETURN
674  end
675 
676 c************************************************************************
677  function dmapa(f, a, b, c, d)
678 c************************************************************************
679 c
680 c.....Derivative of linear map of f from range (a,b) to range (c,d)
681 c.....(silly, but made for some generality purposes...)
682 c
683  implicit none
684 c.....variables
685  double precision a, b, c, d, f, dmapa
686 
687 c.....map vector (no check of bounds...)
688  dmapa = ( d - c ) / ( b - a )
689 
690  RETURN
691  end
692 
693 c************************************************************************
694  subroutine vmapa(vin, n, a, b, c, d, vout)
695 c************************************************************************
696 c
697 c.....Linear map of vector elements from range (a,b) to range (c,d)
698 c
699  implicit none
700 c.....variables
701  INTEGER n
702  double precision a, b, c, d, vin(n), vout(n)
703 c.....helper variables
704  integer k
705  double precision t
706 c.....f2py statements
707 cf2py intent(out) vout
708 
709 c.....map vector (no check of bounds...)
710  t = ( d - c ) / ( b - a )
711  do k=1,n
712  vout(k) = c + ( vin(k) - a ) * t
713  enddo
714 
715  RETURN
716  end
717 
718 c************************************************************************
719  subroutine mmapa(mmin, m, n, a, b, c, d, mmout)
720 c************************************************************************
721 c
722 c.....Linear map of matrix elements from range (a,b) to range (c,d)
723 c
724  implicit none
725 c.....variables
726  INTEGER m,n
727  double precision a, b, c, d, mmin(m, n), mmout(m,n)
728 c.....helper variables
729  integer j,k
730  double precision t
731 c.....f2py statements
732 cf2py intent(out) mmout
733 
734 c.....map matrix (no check of bounds...)
735  t = ( d - c ) / ( b - a )
736  do j=1,m
737  do k=1,n
738  mmout(j,k) = c + ( mmin(j,k) - a ) * t
739  enddo
740  enddo
741 
742  RETURN
743  end
744 
745 
subroutine grad(x, conec, n, bconecno, bn, units, u, inno, i, outno, o, Input, Targ, p, xprime)
Definition: ffnet.f:103
subroutine vmapa(vin, n, a, b, c, d, vout)
Definition: ffnet.f:695
double precision function mapa(f, a, b, c, d)
Definition: ffnet.f:659
subroutine normcall2(x, conec, n, units, u, inno, i, outno, o, eni, deo, input, p, output)
Definition: ffnet.f:445
#define sign(x)
Definition: misc.h:95
subroutine sqerror(x, conec, n, units, u, inno, i, outno, o, Input, Targ, p, sqerr)
Definition: ffnet.f:68
double precision function dmapa(f, a, b, c, d)
Definition: ffnet.f:678
subroutine getout(units, u, outno, o, deo, output)
Definition: ffnet.f:637
subroutine pikaiaff(x, ffn, conec, n, units, u, inno, i, outno, o, Input, Targ, p, bound1, bound2, isqerr)
Definition: ffnet.f:310
subroutine momentum(x, conec, n, bconecno, bn, units, u, inno, i, outno, o, Input, Targ, p, eta, moment, maxiter)
Definition: ffnet.f:522
subroutine setin(input, inno, i, eni, units, u)
Definition: ffnet.f:615
subroutine recall(x, conec, n, units, u, inno, i, outno, o, input, output)
Definition: ffnet.f:175
subroutine func(x, conec, n, bconecno, bn, units, u, inno, i, outno, o, Input, Targ, p, sqerr)
Definition: ffnet.f:287
#define max(A, B)
Definition: main_biosmap.c:61
subroutine normdiff2(x, conec, n, dconecno, dn, dconecmk, units, u, inno, i, outno, o, eni, ded, input, p, deriv)
Definition: ffnet.f:478
subroutine diff(x, conec, n, dconecno, dn, dconecmk, units, u, inno, i, outno, o, input, deriv)
Definition: ffnet.f:205
#define min(A, B)
Definition: main_biosmap.c:62
subroutine normcall(x, conec, n, units, u, inno, i, outno, o, eni, deo, input, output)
Definition: ffnet.f:339
subroutine rprop(x, conec, n, bconecno, bn, units, u, inno, i, outno, o, Input, Targ, p, a, b, mimin, mimax, xmi, maxiter)
Definition: ffnet.f:561
subroutine normdiff(x, conec, n, dconecno, dn, dconecmk, units, u, inno, i, outno, o, eni, ded, input, deriv)
Definition: ffnet.f:367
subroutine prop(x, conec, n, units, u)
Definition: ffnet.f:19
#define f
Definition: l1_czcs_hdf.c:702
subroutine mmapa(mmin, m, n, a, b, c, d, mmout)
Definition: ffnet.f:720