NEURON
nvector_nrnserial_ld.cpp
Go to the documentation of this file.
1 /*
2  * -----------------------------------------------------------------
3  * $Revision: 2212 $
4  * $Date: 2008-09-08 16:32:26 +0200 (Mon, 08 Sep 2008) $
5  * -----------------------------------------------------------------
6  * Programmer(s): Scott D. Cohen, Alan C. Hindmarsh, Radu Serban,
7  * and Aaron Collier @ LLNL
8  * -----------------------------------------------------------------
9  * Copyright (c) 2002, The Regents of the University of California.
10  * Produced at the Lawrence Livermore National Laboratory.
11  * All rights reserved.
12  * For details, see sundials/shared/LICENSE.
13  * -----------------------------------------------------------------
14  * This is the implementation file for a serial implementation
15  * of the NVECTOR package.
16  * -----------------------------------------------------------------
17  */
18 
19 #define USELONGDOUBLE 1
20 
21 #include <../../nrnconf.h>
22 #include <hocassrt.h>
23 #include <nrnassrt.h>
24 #if HAVE_POSIX_MEMALIGN
25 #define HAVE_MEMALIGN 1
26 #endif
27 #if HAVE_MEMALIGN
28 #undef _XOPEN_SOURCE /* avoid warnings about redefining this */
29 #define _XOPEN_SOURCE 600
30 #endif
31 
32 #include <stdio.h>
33 #include <stdlib.h>
34 
35 #include "nvector_nrnserial_ld.h"
36 #include "shared/sundialsmath.h"
37 #include "shared/sundialstypes.h"
38 
39 #define ZERO RCONST(0.0)
40 #define HALF RCONST(0.5)
41 #define ONE RCONST(1.0)
42 #define ONEPT5 RCONST(1.5)
43 
44 #if USELONGDOUBLE
45 #define ldrealtype long double
46 #else
47 #define ldrealtype realtype
48 #endif
49 
50 /* Private function prototypes */
51 /* z=x */
52 static void VCopy_NrnSerialLD(N_Vector x, N_Vector z);
53 /* z=x+y */
54 static void VSum_NrnSerialLD(N_Vector x, N_Vector y, N_Vector z);
55 /* z=x-y */
56 static void VDiff_NrnSerialLD(N_Vector x, N_Vector y, N_Vector z);
57 /* z=-x */
58 static void VNeg_NrnSerialLD(N_Vector x, N_Vector z);
59 /* z=c(x+y) */
60 static void VScaleSum_NrnSerialLD(realtype c, N_Vector x, N_Vector y, N_Vector z);
61 /* z=c(x-y) */
62 static void VScaleDiff_NrnSerialLD(realtype c, N_Vector x, N_Vector y, N_Vector z);
63 /* z=ax+y */
64 static void VLin1_NrnSerialLD(realtype a, N_Vector x, N_Vector y, N_Vector z);
65 /* z=ax-y */
66 static void VLin2_NrnSerialLD(realtype a, N_Vector x, N_Vector y, N_Vector z);
67 /* y <- ax+y */
68 static void Vaxpy_NrnSerialLD(realtype a, N_Vector x, N_Vector y);
69 /* x <- ax */
70 static void VScaleBy_NrnSerialLD(realtype a, N_Vector x);
71 
72 /*
73  * -----------------------------------------------------------------
74  * exported functions
75  * -----------------------------------------------------------------
76  */
77 
78 /* ----------------------------------------------------------------------------
79  * Function to create a new empty serial vector
80  */
81 
82 N_Vector N_VNewEmpty_NrnSerialLD(long int length)
83 {
84  N_Vector v;
85  N_Vector_Ops ops;
87 
88  /* Create vector */
89  v = (N_Vector) malloc(sizeof *v);
90  if (v == NULL) return(NULL);
91 
92  /* Create vector operation structure */
93  ops = (N_Vector_Ops) malloc(sizeof(struct _generic_N_Vector_Ops));
94  if (ops == NULL) {free(v);return(NULL);}
95 
96  ops->nvclone = N_VClone_NrnSerialLD;
97  ops->nvdestroy = N_VDestroy_NrnSerialLD;
98  ops->nvspace = N_VSpace_NrnSerialLD;
99  ops->nvgetarraypointer = N_VGetArrayPointer_NrnSerialLD;
100  ops->nvsetarraypointer = N_VSetArrayPointer_NrnSerialLD;
101  ops->nvlinearsum = N_VLinearSum_NrnSerialLD;
102  ops->nvconst = N_VConst_NrnSerialLD;
103  ops->nvprod = N_VProd_NrnSerialLD;
104  ops->nvdiv = N_VDiv_NrnSerialLD;
105  ops->nvscale = N_VScale_NrnSerialLD;
106  ops->nvabs = N_VAbs_NrnSerialLD;
107  ops->nvinv = N_VInv_NrnSerialLD;
108  ops->nvaddconst = N_VAddConst_NrnSerialLD;
109  ops->nvdotprod = N_VDotProd_NrnSerialLD;
110  ops->nvmaxnorm = N_VMaxNorm_NrnSerialLD;
111  ops->nvwrmsnormmask = N_VWrmsNormMask_NrnSerialLD;
112  ops->nvwrmsnorm = N_VWrmsNorm_NrnSerialLD;
113  ops->nvmin = N_VMin_NrnSerialLD;
114  ops->nvwl2norm = N_VWL2Norm_NrnSerialLD;
115  ops->nvl1norm = N_VL1Norm_NrnSerialLD;
116  ops->nvcompare = N_VCompare_NrnSerialLD;
117  ops->nvinvtest = N_VInvTest_NrnSerialLD;
118  ops->nvconstrmask = N_VConstrMask_NrnSerialLD;
119  ops->nvminquotient = N_VMinQuotient_NrnSerialLD;
120 
121  /* Create content */
122  content = (N_VectorContent_NrnSerialLD) malloc(sizeof(struct _N_VectorContent_NrnSerialLD));
123  if (content == NULL) {free(ops);free(v);return(NULL);}
124 
125  content->length = length;
126  content->own_data = FALSE;
127  content->data = NULL;
128 
129  /* Attach content and ops */
130  v->content = content;
131  v->ops = ops;
132 
133  return(v);
134 }
135 
136 /* ----------------------------------------------------------------------------
137  * Function to create a new serial vector
138  */
139 
140 N_Vector N_VNew_NrnSerialLD(long int length)
141 {
142  N_Vector v;
143  realtype *data;
144 
145  v = N_VNewEmpty_NrnSerialLD(length);
146  if (v == NULL) return(NULL);
147 
148  /* Create data */
149  if (length > 0) {
150 
151  /* Allocate memory */
152 #if HAVE_MEMALIGN
153  nrn_assert(posix_memalign((void**)&data, 64, length*sizeof(realtype)) == 0);
154 #else
155  data = (realtype *) malloc(length * sizeof(realtype));
156 #endif
157  if(data == NULL) {N_VDestroy_NrnSerialLD(v);return(NULL);}
158 
159  /* Attach data */
160  NV_OWN_DATA_S_LD(v) = TRUE;
161  NV_DATA_S_LD(v) = data;
162 
163  }
164 
165  return(v);
166 }
167 
168 /* ----------------------------------------------------------------------------
169  * Function to clone from a template a new vector with empty (NULL) data array
170  */
171 
172 N_Vector N_VCloneEmpty_NrnSerialLD(N_Vector w)
173 {
174  N_Vector v;
175  N_Vector_Ops ops;
177 
178  if (w == NULL) return(NULL);
179 
180  /* Create vector */
181  v = (N_Vector) malloc(sizeof *v);
182  if (v == NULL) return(NULL);
183 
184  /* Create vector operation structure */
185  ops = (N_Vector_Ops) malloc(sizeof(struct _generic_N_Vector_Ops));
186  if (ops == NULL) {free(v);return(NULL);}
187 
188  ops->nvclone = w->ops->nvclone;
189  ops->nvdestroy = w->ops->nvdestroy;
190  ops->nvspace = w->ops->nvspace;
191  ops->nvgetarraypointer = w->ops->nvgetarraypointer;
192  ops->nvsetarraypointer = w->ops->nvsetarraypointer;
193  ops->nvlinearsum = w->ops->nvlinearsum;
194  ops->nvconst = w->ops->nvconst;
195  ops->nvprod = w->ops->nvprod;
196  ops->nvdiv = w->ops->nvdiv;
197  ops->nvscale = w->ops->nvscale;
198  ops->nvabs = w->ops->nvabs;
199  ops->nvinv = w->ops->nvinv;
200  ops->nvaddconst = w->ops->nvaddconst;
201  ops->nvdotprod = w->ops->nvdotprod;
202  ops->nvmaxnorm = w->ops->nvmaxnorm;
203  ops->nvwrmsnormmask = w->ops->nvwrmsnormmask;
204  ops->nvwrmsnorm = w->ops->nvwrmsnorm;
205  ops->nvmin = w->ops->nvmin;
206  ops->nvwl2norm = w->ops->nvwl2norm;
207  ops->nvl1norm = w->ops->nvl1norm;
208  ops->nvcompare = w->ops->nvcompare;
209  ops->nvinvtest = w->ops->nvinvtest;
210  ops->nvconstrmask = w->ops->nvconstrmask;
211  ops->nvminquotient = w->ops->nvminquotient;
212 
213  /* Create content */
214  content = (N_VectorContent_NrnSerialLD) malloc(sizeof(struct _N_VectorContent_NrnSerialLD));
215  if (content == NULL) {free(ops);free(v);return(NULL);}
216 
217  content->length = NV_LENGTH_S_LD(w);
218  content->own_data = FALSE;
219  content->data = NULL;
220 
221  /* Attach content and ops */
222  v->content = content;
223  v->ops = ops;
224 
225  return(v);
226 }
227 
228 /* ----------------------------------------------------------------------------
229  * Function to create a serial N_Vector with user data component
230  */
231 
232 N_Vector N_VMake_NrnSerialLD(long int length, realtype *v_data)
233 {
234  N_Vector v;
235 
236  v = N_VNewEmpty_NrnSerialLD(length);
237  if (v == NULL) return(NULL);
238 
239  if (length > 0) {
240  /* Attach data */
241  NV_OWN_DATA_S_LD(v) = FALSE;
242  NV_DATA_S_LD(v) = v_data;
243  }
244 
245  return(v);
246 }
247 
248 /* ----------------------------------------------------------------------------
249  * Function to create an array of new serial vectors.
250  */
251 
252 N_Vector *N_VNewVectorArray_NrnSerialLD(int count, long int length)
253 {
254  N_Vector *vs;
255  int j;
256 
257  if (count <= 0) return(NULL);
258 
259  vs = (N_Vector *) malloc(count * sizeof(N_Vector));
260  if(vs == NULL) return(NULL);
261 
262  for (j=0; j<count; j++) {
263  vs[j] = N_VNew_NrnSerialLD(length);
264  if (vs[j] == NULL) {
266  return(NULL);
267  }
268  }
269 
270  return(vs);
271 }
272 
273 /* ----------------------------------------------------------------------------
274  * Function to create an array of new serial vectors with NULL data array.
275  */
276 
277 N_Vector *N_VNewVectorArrayEmpty_NrnSerialLD(int count, long int length)
278 {
279  N_Vector *vs;
280  int j;
281 
282  if (count <= 0) return(NULL);
283 
284  vs = (N_Vector *) malloc(count * sizeof(N_Vector));
285  if(vs == NULL) return(NULL);
286 
287  for (j=0; j<count; j++) {
288  vs[j] = N_VNewEmpty_NrnSerialLD(length);
289  if (vs[j] == NULL) {
291  return(NULL);
292  }
293  }
294 
295  return(vs);
296 }
297 
298 /* ----------------------------------------------------------------------------
299  * Function to free an array created with N_VNewVectorArray_NrnSerialLD
300  */
301 
302 void N_VDestroyVectorArray_NrnSerialLD(N_Vector *vs, int count)
303 {
304  int j;
305 
306  for (j = 0; j < count; j++) N_VDestroy_NrnSerialLD(vs[j]);
307 
308  free(vs);
309 }
310 
311 /* ----------------------------------------------------------------------------
312  * Function to print the a serial vector
313  */
314 
315 void N_VPrint_NrnSerialLD(N_Vector x)
316 {
317  long int i, N;
318  realtype *xd;
319 
320  N = NV_LENGTH_S_LD(x);
321  xd = NV_DATA_S_LD(x);
322 
323  for (i=0; i < N; i++) {
324 #if defined(SUNDIALS_EXTENDED_PRECISION)
325  printf("%11.8Lg\n", *xd++);
326 #elif defined(SUNDIALS_DOUBLE_PRECISION)
327  printf("%11.8lg\n", *xd++);
328 #else
329  printf("%11.8g\n", *xd++);
330 #endif
331  }
332  printf("\n");
333 }
334 
335 /*
336  * -----------------------------------------------------------------
337  * implementation of vector operations
338  * -----------------------------------------------------------------
339  */
340 
341 N_Vector N_VClone_NrnSerialLD(N_Vector w)
342 {
343  N_Vector v;
344  realtype *data;
345  long int length;
346 
348  if (v == NULL) return(NULL);
349 
350  length = NV_LENGTH_S_LD(w);
351 
352  /* Create data */
353  if (length > 0) {
354 
355  /* Allocate memory */
356 #if HAVE_MEMALIGN
357  nrn_assert(posix_memalign((void**)&data, 64, length*sizeof(realtype)) == 0);
358 #else
359  data = (realtype *) malloc(length * sizeof(realtype));
360 #endif
361  if(data == NULL) {N_VDestroy_NrnSerialLD(v);return(NULL);}
362 
363  /* Attach data */
364  NV_OWN_DATA_S_LD(v) = TRUE;
365  NV_DATA_S_LD(v) = data;
366 
367  }
368 
369  return(v);
370 }
371 
372 void N_VDestroy_NrnSerialLD(N_Vector v)
373 {
374  if (NV_OWN_DATA_S_LD(v) == TRUE)
375  free(NV_DATA_S_LD(v));
376  free(v->content);
377  free(v->ops);
378  free(v);
379 }
380 
381 void N_VSpace_NrnSerialLD(N_Vector v, long int *lrw, long int *liw)
382 {
383  *lrw = NV_LENGTH_S_LD(v);
384  *liw = 1;
385 }
386 
388 {
389  realtype *v_data;
390 
391  v_data = NV_DATA_S_LD(v);
392 
393  return(v_data);
394 }
395 
396 void N_VSetArrayPointer_NrnSerialLD(realtype *v_data, N_Vector v)
397 {
398  if (NV_LENGTH_S_LD(v) > 0) NV_DATA_S_LD(v) = v_data;
399 }
400 
401 void N_VLinearSum_NrnSerialLD(realtype a, N_Vector x, realtype b, N_Vector y, N_Vector z)
402 {
403  long int i, N;
404  realtype c, *xd, *yd, *zd;
405  N_Vector v1, v2;
406  booleantype test;
407 
408  if ((b == ONE) && (z == y)) { /* BLAS usage: axpy y <- ax+y */
409  Vaxpy_NrnSerialLD(a,x,y);
410  return;
411  }
412 
413  if ((a == ONE) && (z == x)) { /* BLAS usage: axpy x <- by+x */
414  Vaxpy_NrnSerialLD(b,y,x);
415  return;
416  }
417 
418  /* Case: a == b == 1.0 */
419 
420  if ((a == ONE) && (b == ONE)) {
421  VSum_NrnSerialLD(x, y, z);
422  return;
423  }
424 
425  /* Cases: (1) a == 1.0, b = -1.0, (2) a == -1.0, b == 1.0 */
426 
427  if ((test = ((a == ONE) && (b == -ONE))) || ((a == -ONE) && (b == ONE))) {
428  v1 = test ? y : x;
429  v2 = test ? x : y;
430  VDiff_NrnSerialLD(v2, v1, z);
431  return;
432  }
433 
434  /* Cases: (1) a == 1.0, b == other or 0.0, (2) a == other or 0.0, b == 1.0 */
435  /* if a or b is 0.0, then user should have called N_VScale */
436 
437  if ((test = (a == ONE)) || (b == ONE)) {
438  c = test ? b : a;
439  v1 = test ? y : x;
440  v2 = test ? x : y;
441  VLin1_NrnSerialLD(c, v1, v2, z);
442  return;
443  }
444 
445  /* Cases: (1) a == -1.0, b != 1.0, (2) a != 1.0, b == -1.0 */
446 
447  if ((test = (a == -ONE)) || (b == -ONE)) {
448  c = test ? b : a;
449  v1 = test ? y : x;
450  v2 = test ? x : y;
451  VLin2_NrnSerialLD(c, v1, v2, z);
452  return;
453  }
454 
455  /* Case: a == b */
456  /* catches case both a and b are 0.0 - user should have called N_VConst */
457 
458  if (a == b) {
459  VScaleSum_NrnSerialLD(a, x, y, z);
460  return;
461  }
462 
463  /* Case: a == -b */
464 
465  if (a == -b) {
466  VScaleDiff_NrnSerialLD(a, x, y, z);
467  return;
468  }
469 
470  /* Do all cases not handled above:
471  (1) a == other, b == 0.0 - user should have called N_VScale
472  (2) a == 0.0, b == other - user should have called N_VScale
473  (3) a,b == other, a !=b, a != -b */
474 
475  N = NV_LENGTH_S_LD(x);
476  xd = NV_DATA_S_LD(x);
477  yd = NV_DATA_S_LD(y);
478  zd = NV_DATA_S_LD(z);
479 
480  for (i=0; i < N; i++)
481  *zd++ = a * (*xd++) + b * (*yd++);
482 }
483 
484 void N_VConst_NrnSerialLD(realtype c, N_Vector z)
485 {
486  long int i, N;
487  realtype *zd;
488 
489  N = NV_LENGTH_S_LD(z);
490  zd = NV_DATA_S_LD(z);
491 
492  for (i=0; i < N; i++)
493  *zd++ = c;
494 }
495 
496 void N_VProd_NrnSerialLD(N_Vector x, N_Vector y, N_Vector z)
497 {
498  long int i, N;
499  realtype *xd, *yd, *zd;
500 
501  N = NV_LENGTH_S_LD(x);
502  xd = NV_DATA_S_LD(x);
503  yd = NV_DATA_S_LD(y);
504  zd = NV_DATA_S_LD(z);
505 
506  for (i=0; i < N; i++)
507  *zd++ = (*xd++) * (*yd++);
508 }
509 
510 void N_VDiv_NrnSerialLD(N_Vector x, N_Vector y, N_Vector z)
511 {
512  long int i, N;
513  realtype *xd, *yd, *zd;
514 
515  N = NV_LENGTH_S_LD(x);
516  xd = NV_DATA_S_LD(x);
517  yd = NV_DATA_S_LD(y);
518  zd = NV_DATA_S_LD(z);
519 
520  for (i=0; i < N; i++)
521  *zd++ = (*xd++) / (*yd++);
522 }
523 
524 void N_VScale_NrnSerialLD(realtype c, N_Vector x, N_Vector z)
525 {
526  long int i, N;
527  realtype *xd, *zd;
528 
529  if (z == x) { /* BLAS usage: scale x <- cx */
530  VScaleBy_NrnSerialLD(c, x);
531  return;
532  }
533 
534  if (c == ONE) {
535  VCopy_NrnSerialLD(x, z);
536  } else if (c == -ONE) {
537  VNeg_NrnSerialLD(x, z);
538  } else {
539  N = NV_LENGTH_S_LD(x);
540  xd = NV_DATA_S_LD(x);
541  zd = NV_DATA_S_LD(z);
542  for (i=0; i < N; i++)
543  *zd++ = c * (*xd++);
544  }
545 }
546 
547 void N_VAbs_NrnSerialLD(N_Vector x, N_Vector z)
548 {
549  long int i, N;
550  realtype *xd, *zd;
551 
552  N = NV_LENGTH_S_LD(x);
553  xd = NV_DATA_S_LD(x);
554  zd = NV_DATA_S_LD(z);
555 
556  for (i=0; i < N; i++, xd++, zd++)
557  *zd = ABS(*xd);
558 }
559 
560 void N_VInv_NrnSerialLD(N_Vector x, N_Vector z)
561 {
562  long int i, N;
563  realtype *xd, *zd;
564 
565  N = NV_LENGTH_S_LD(x);
566  xd = NV_DATA_S_LD(x);
567  zd = NV_DATA_S_LD(z);
568 
569  for (i=0; i < N; i++)
570  *zd++ = ONE / (*xd++);
571 }
572 
573 void N_VAddConst_NrnSerialLD(N_Vector x, realtype b, N_Vector z)
574 {
575  long int i, N;
576  realtype *xd, *zd;
577 
578  N = NV_LENGTH_S_LD(x);
579  xd = NV_DATA_S_LD(x);
580  zd = NV_DATA_S_LD(z);
581 
582  for (i=0; i < N; i++)
583  *zd++ = (*xd++) + b;
584 }
585 
586 realtype N_VDotProd_NrnSerialLD(N_Vector x, N_Vector y)
587 {
588  long int i, N;
589  realtype sum = ZERO, *xd, *yd;
590 
591  N = NV_LENGTH_S_LD(x);
592  xd = NV_DATA_S_LD(x);
593  yd = NV_DATA_S_LD(y);
594 
595  for (i=0; i < N; i++)
596  sum += (*xd++) * (*yd++);
597 
598  return(sum);
599 }
600 
601 realtype N_VMaxNorm_NrnSerialLD(N_Vector x)
602 {
603  long int i, N;
604  realtype max = ZERO, *xd;
605 
606  N = NV_LENGTH_S_LD(x);
607  xd = NV_DATA_S_LD(x);
608 
609  for (i=0; i < N; i++, xd++) {
610  if (ABS(*xd) > max) max = ABS(*xd);
611  }
612 
613  return(max);
614 }
615 
616 realtype N_VWrmsNorm_NrnSerialLD(N_Vector x, N_Vector w)
617 {
618  long int i, N;
619  realtype prodi, *xd, *wd;
620  ldrealtype sum = ZERO;
621 
622  N = NV_LENGTH_S_LD(x);
623  xd = NV_DATA_S_LD(x);
624  wd = NV_DATA_S_LD(w);
625 
626  for (i=0; i < N; i++) {
627  prodi = (*xd++) * (*wd++);
628  sum += prodi * prodi;
629  }
630 
631  return(RSqrt((realtype)sum / N));
632 }
633 
634 realtype N_VWrmsNormMask_NrnSerialLD(N_Vector x, N_Vector w, N_Vector id)
635 {
636  long int i, N;
637  realtype prodi, *xd, *wd, *idd;
638  ldrealtype sum = ZERO;
639 
640  N = NV_LENGTH_S_LD(x);
641  xd = NV_DATA_S_LD(x);
642  wd = NV_DATA_S_LD(w);
643  idd = NV_DATA_S_LD(id);
644 
645  for (i=0; i < N; i++) {
646  if (idd[i] > ZERO) {
647  prodi = xd[i] * wd[i];
648  sum += prodi * prodi;
649  }
650  }
651 
652  return(RSqrt((realtype)sum / N));
653 }
654 
655 realtype N_VMin_NrnSerialLD(N_Vector x)
656 {
657  long int i, N;
658  realtype min, *xd;
659 
660  N = NV_LENGTH_S_LD(x);
661  xd = NV_DATA_S_LD(x);
662 
663  min = xd[0];
664 
665  xd++;
666  for (i=1; i < N; i++, xd++) {
667  if ((*xd) < min) min = *xd;
668  }
669 
670  return(min);
671 }
672 
673 realtype N_VWL2Norm_NrnSerialLD(N_Vector x, N_Vector w)
674 {
675  long int i, N;
676  realtype prodi, *xd, *wd;
677  ldrealtype sum = ZERO;
678 
679  N = NV_LENGTH_S_LD(x);
680  xd = NV_DATA_S_LD(x);
681  wd = NV_DATA_S_LD(w);
682 
683  for (i=0; i < N; i++) {
684  prodi = (*xd++) * (*wd++);
685  sum += prodi * prodi;
686  }
687 
688  return(RSqrt((realtype)sum));
689 }
690 
691 realtype N_VL1Norm_NrnSerialLD(N_Vector x)
692 {
693  long int i, N;
694  realtype *xd;
695  ldrealtype sum = ZERO;
696 
697  N = NV_LENGTH_S_LD(x);
698  xd = NV_DATA_S_LD(x);
699 
700  for (i=0; i<N; i++)
701  sum += ABS(xd[i]);
702 
703  return((realtype)sum);
704 }
705 
706 void N_VOneMask_NrnSerialLD(N_Vector x)
707 {
708  long int i, N;
709  realtype *xd;
710 
711  N = NV_LENGTH_S_LD(x);
712  xd = NV_DATA_S_LD(x);
713 
714  for (i=0; i<N; i++,xd++) {
715  if (*xd != ZERO) *xd = ONE;
716  }
717 }
718 
719 void N_VCompare_NrnSerialLD(realtype c, N_Vector x, N_Vector z)
720 {
721  long int i, N;
722  realtype *xd, *zd;
723 
724  N = NV_LENGTH_S_LD(x);
725  xd = NV_DATA_S_LD(x);
726  zd = NV_DATA_S_LD(z);
727 
728  for (i=0; i < N; i++, xd++, zd++) {
729  *zd = (ABS(*xd) >= c) ? ONE : ZERO;
730  }
731 }
732 
733 booleantype N_VInvTest_NrnSerialLD(N_Vector x, N_Vector z)
734 {
735  long int i, N;
736  realtype *xd, *zd;
737 
738  N = NV_LENGTH_S_LD(x);
739  xd = NV_DATA_S_LD(x);
740  zd = NV_DATA_S_LD(z);
741 
742  for (i=0; i < N; i++) {
743  if (*xd == ZERO) return(FALSE);
744  *zd++ = ONE / (*xd++);
745  }
746 
747  return(TRUE);
748 }
749 
750 booleantype N_VConstrMask_NrnSerialLD(N_Vector c, N_Vector x, N_Vector m)
751 {
752  long int i, N;
753  booleantype test;
754  realtype *cd, *xd, *md;
755 
756  N = NV_LENGTH_S_LD(x);
757  xd = NV_DATA_S_LD(x);
758  cd = NV_DATA_S_LD(c);
759  md = NV_DATA_S_LD(m);
760 
761  test = TRUE;
762 
763  for (i=0; i<N; i++, cd++, xd++, md++) {
764  *md = ZERO;
765  if (*cd == ZERO) continue;
766  if (*cd > ONEPT5 || (*cd) < -ONEPT5) {
767  if ( (*xd)*(*cd) <= ZERO) {test = FALSE; *md = ONE; }
768  continue;
769  }
770  if ( (*cd) > HALF || (*cd) < -HALF) {
771  if ( (*xd)*(*cd) < ZERO ) {test = FALSE; *md = ONE; }
772  }
773  }
774  return(test);
775 }
776 
777 realtype N_VMinQuotient_NrnSerialLD(N_Vector num, N_Vector denom)
778 {
779  booleantype notEvenOnce;
780  long int i, N;
781  realtype *nd, *dd, min=0.0;
782 
783  N = NV_LENGTH_S_LD(num);
784  nd = NV_DATA_S_LD(num);
785  dd = NV_DATA_S_LD(denom);
786 
787  notEvenOnce = TRUE;
788 
789  for (i = 0; i < N; i++, nd++, dd++) {
790  if (*dd == ZERO) continue;
791  else {
792  if (notEvenOnce) {
793  min = *nd / *dd ;
794  notEvenOnce = FALSE;
795  }
796  else min = MIN(min, (*nd) / (*dd));
797  }
798  }
799 
800  if (notEvenOnce || (N == 0)) min = BIG_REAL;
801 
802  return(min);
803 }
804 
805 /*
806  * -----------------------------------------------------------------
807  * private functions
808  * -----------------------------------------------------------------
809  */
810 
811 static void VCopy_NrnSerialLD(N_Vector x, N_Vector z)
812 {
813  long int i, N;
814  realtype *xd, *zd;
815 
816  N = NV_LENGTH_S_LD(x);
817  xd = NV_DATA_S_LD(x);
818  zd = NV_DATA_S_LD(z);
819 
820  for (i=0; i < N; i++)
821  *zd++ = *xd++;
822 }
823 
824 static void VSum_NrnSerialLD(N_Vector x, N_Vector y, N_Vector z)
825 {
826  long int i, N;
827  realtype *xd, *yd, *zd;
828 
829  N = NV_LENGTH_S_LD(x);
830  xd = NV_DATA_S_LD(x);
831  yd = NV_DATA_S_LD(y);
832  zd = NV_DATA_S_LD(z);
833 
834  for (i=0; i < N; i++)
835  *zd++ = (*xd++) + (*yd++);
836 }
837 
838 static void VDiff_NrnSerialLD(N_Vector x, N_Vector y, N_Vector z)
839 {
840  long int i, N;
841  realtype *xd, *yd, *zd;
842 
843  N = NV_LENGTH_S_LD(x);
844  xd = NV_DATA_S_LD(x);
845  yd = NV_DATA_S_LD(y);
846  zd = NV_DATA_S_LD(z);
847 
848  for (i=0; i < N; i++)
849  *zd++ = (*xd++) - (*yd++);
850 }
851 
852 static void VNeg_NrnSerialLD(N_Vector x, N_Vector z)
853 {
854  long int i, N;
855  realtype *xd, *zd;
856 
857  N = NV_LENGTH_S_LD(x);
858  xd = NV_DATA_S_LD(x);
859  zd = NV_DATA_S_LD(z);
860 
861  for (i=0; i < N; i++)
862  *zd++ = -(*xd++);
863 }
864 
865 static void VScaleSum_NrnSerialLD(realtype c, N_Vector x, N_Vector y, N_Vector z)
866 {
867  long int i, N;
868  realtype *xd, *yd, *zd;
869 
870  N = NV_LENGTH_S_LD(x);
871  xd = NV_DATA_S_LD(x);
872  yd = NV_DATA_S_LD(y);
873  zd = NV_DATA_S_LD(z);
874 
875  for (i=0; i < N; i++)
876  *zd++ = c * ((*xd++) + (*yd++));
877 }
878 
879 static void VScaleDiff_NrnSerialLD(realtype c, N_Vector x, N_Vector y, N_Vector z)
880 {
881  long int i, N;
882  realtype *xd, *yd, *zd;
883 
884  N = NV_LENGTH_S_LD(x);
885  xd = NV_DATA_S_LD(x);
886  yd = NV_DATA_S_LD(y);
887  zd = NV_DATA_S_LD(z);
888 
889  for (i=0; i < N; i++)
890  *zd++ = c * ((*xd++) - (*yd++));
891 }
892 
893 static void VLin1_NrnSerialLD(realtype a, N_Vector x, N_Vector y, N_Vector z)
894 {
895  long int i, N;
896  realtype *xd, *yd, *zd;
897 
898  N = NV_LENGTH_S_LD(x);
899  xd = NV_DATA_S_LD(x);
900  yd = NV_DATA_S_LD(y);
901  zd = NV_DATA_S_LD(z);
902 
903  for (i=0; i < N; i++)
904  *zd++ = a * (*xd++) + (*yd++);
905 }
906 
907 static void VLin2_NrnSerialLD(realtype a, N_Vector x, N_Vector y, N_Vector z)
908 {
909  long int i, N;
910  realtype *xd, *yd, *zd;
911 
912  N = NV_LENGTH_S_LD(x);
913  xd = NV_DATA_S_LD(x);
914  yd = NV_DATA_S_LD(y);
915  zd = NV_DATA_S_LD(z);
916 
917  for (i=0; i < N; i++)
918  *zd++ = a * (*xd++) - (*yd++);
919 }
920 
921 static void Vaxpy_NrnSerialLD(realtype a, N_Vector x, N_Vector y)
922 {
923  long int i, N;
924  realtype *xd, *yd;
925 
926  N = NV_LENGTH_S_LD(x);
927  xd = NV_DATA_S_LD(x);
928  yd = NV_DATA_S_LD(y);
929 
930  if (a == ONE) {
931  for (i=0; i < N; i++)
932  *yd++ += (*xd++);
933  return;
934  }
935 
936  if (a == -ONE) {
937  for (i=0; i < N; i++)
938  *yd++ -= (*xd++);
939  return;
940  }
941 
942  for (i=0; i < N; i++)
943  *yd++ += a * (*xd++);
944 }
945 
946 static void VScaleBy_NrnSerialLD(realtype a, N_Vector x)
947 {
948  long int i, N;
949  realtype *xd;
950 
951  N = NV_LENGTH_S_LD(x);
952  xd = NV_DATA_S_LD(x);
953 
954  for (i=0; i < N; i++)
955  *xd++ *= a;
956 }
#define data
Definition: rbtqueue.cpp:49
#define nrn_assert(ex)
Definition: nrnassrt.h:35
double max(double a, double b)
Definition: geometry3d.cpp:22
static void VLin1_NrnSerialLD(realtype a, N_Vector x, N_Vector y, N_Vector z)
static void VNeg_NrnSerialLD(N_Vector x, N_Vector z)
void N_VDiv_NrnSerialLD(N_Vector x, N_Vector y, N_Vector z)
N_Vector N_VClone_NrnSerialLD(N_Vector w)
static void VCopy_NrnSerialLD(N_Vector x, N_Vector z)
#define ONE
static void VScaleDiff_NrnSerialLD(realtype c, N_Vector x, N_Vector y, N_Vector z)
#define min(a, b)
Definition: matrix.h:157
void N_VCompare_NrnSerialLD(realtype c, N_Vector x, N_Vector z)
#define TRUE
Definition: err.c:57
#define NV_LENGTH_S_LD(v)
#define ZERO
void N_VAbs_NrnSerialLD(N_Vector x, N_Vector z)
N_Vector * N_VNewVectorArray_NrnSerialLD(int count, long int length)
nd
Definition: treeset.cpp:893
#define v
Definition: md1redef.h:4
void N_VProd_NrnSerialLD(N_Vector x, N_Vector y, N_Vector z)
realtype * N_VGetArrayPointer_NrnSerialLD(N_Vector v)
N_Vector * N_VNewVectorArrayEmpty_NrnSerialLD(int count, long int length)
realtype N_VL1Norm_NrnSerialLD(N_Vector x)
#define ABS(x)
Definition: extras.c:201
N_Vector N_VCloneEmpty_NrnSerialLD(N_Vector w)
#define HALF
static void VLin2_NrnSerialLD(realtype a, N_Vector x, N_Vector y, N_Vector z)
realtype N_VWL2Norm_NrnSerialLD(N_Vector x, N_Vector w)
static void VScaleSum_NrnSerialLD(realtype c, N_Vector x, N_Vector y, N_Vector z)
static void VSum_NrnSerialLD(N_Vector x, N_Vector y, N_Vector z)
#define ldrealtype
realtype N_VMin_NrnSerialLD(N_Vector x)
realtype N_VMaxNorm_NrnSerialLD(N_Vector x)
struct _N_VectorContent_NrnSerialLD * N_VectorContent_NrnSerialLD
void N_VAddConst_NrnSerialLD(N_Vector x, realtype b, N_Vector z)
booleantype N_VInvTest_NrnSerialLD(N_Vector x, N_Vector z)
#define printf
Definition: mwprefix.h:26
N_Vector N_VMake_NrnSerialLD(long int length, realtype *v_data)
void N_VSpace_NrnSerialLD(N_Vector v, long int *lrw, long int *liw)
N_Vector N_VNew_NrnSerialLD(long int length)
void N_VInv_NrnSerialLD(N_Vector x, N_Vector z)
#define NV_DATA_S_LD(v)
realtype N_VDotProd_NrnSerialLD(N_Vector x, N_Vector y)
size_t j
void N_VScale_NrnSerialLD(realtype c, N_Vector x, N_Vector z)
realtype N_VWrmsNormMask_NrnSerialLD(N_Vector x, N_Vector w, N_Vector id)
void N_VOneMask_NrnSerialLD(N_Vector x)
void N_VPrint_NrnSerialLD(N_Vector x)
#define MIN(a, b)
Definition: grids.h:32
static void VScaleBy_NrnSerialLD(realtype a, N_Vector x)
#define NV_OWN_DATA_S_LD(v)
void N_VDestroy_NrnSerialLD(N_Vector v)
N_Vector N_VNewEmpty_NrnSerialLD(long int length)
realtype N_VWrmsNorm_NrnSerialLD(N_Vector x, N_Vector w)
#define FALSE
Definition: err.c:56
#define i
Definition: md1redef.h:12
#define c
#define ONEPT5
void N_VSetArrayPointer_NrnSerialLD(realtype *v_data, N_Vector v)
void N_VLinearSum_NrnSerialLD(realtype a, N_Vector x, realtype b, N_Vector y, N_Vector z)
void N_VConst_NrnSerialLD(realtype c, N_Vector z)
realtype N_VMinQuotient_NrnSerialLD(N_Vector num, N_Vector denom)
static void VDiff_NrnSerialLD(N_Vector x, N_Vector y, N_Vector z)
return NULL
Definition: cabcode.cpp:461
booleantype N_VConstrMask_NrnSerialLD(N_Vector c, N_Vector x, N_Vector m)
static void Vaxpy_NrnSerialLD(realtype a, N_Vector x, N_Vector y)
void N_VDestroyVectorArray_NrnSerialLD(N_Vector *vs, int count)