NEURON
nvector_nrnparallel_ld.cpp
Go to the documentation of this file.
1 /*
2  * -----------------------------------------------------------------
3  * $Revision: 855 $
4  * $Date: 2005-02-10 00:15:46 +0100 (Thu, 10 Feb 2005) $
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 parallel MPI implementation
15  * of the NVECTOR package.
16  * -----------------------------------------------------------------
17  */
18 
19 #define USELONGDOUBLE 1
20 
21 #include <../../nrnconf.h>
22 #include <hocassrt.h>
23 
24 #include <stdio.h>
25 #include <stdlib.h>
26 
27 /* for NRNMPI_DYNAMICLOAD */
28 #include <nrnmpiuse.h>
29 #if NRNMPI_DYNAMICLOAD
30 extern "C" void nrnmpi_dbl_allreduce_vec(double* src, double* dest, int cnt, int type);
31 extern "C" void nrnmpi_longdbl_allreduce_vec(long double* src,
32  long double* dest,
33  int cnt,
34  int type);
35 extern "C" void nrnmpi_long_allreduce_vec(long* src, long* dest, int cnt, int type);
36 extern int nrnmpi_numprocs;
37 #endif
38 
39 #include "nvector_nrnparallel_ld.h"
40 #if NRNMPI_DYNAMICLOAD
41 #else
42 extern MPI_Comm nrnmpi_comm;
43 #endif
44 #include "sundialsmath.h"
45 #include "sundialstypes.h"
46 
47 #define ZERO RCONST(0.0)
48 #define HALF RCONST(0.5)
49 #define ONE RCONST(1.0)
50 #define ONEPT5 RCONST(1.5)
51 
52 #if USELONGDOUBLE
53 #define ldrealtype long double
54 #else
55 #define ldrealtype realtype
56 #endif
57 
58 
59 /* Error Message */
60 
61 #define BAD_N1 "N_VNew_NrnParallelLD -- Sum of local vector lengths differs from "
62 #define BAD_N2 "input global length. \n\n"
63 #define BAD_N BAD_N1 BAD_N2
64 
65 /* Private function prototypes */
66 
67 /* Reduction operations add/max/min over the processor group */
68 static realtype VAllReduce_NrnParallelLD(realtype d, int op, MPI_Comm comm);
69 #if USELONGDOUBLE
70 /* only for sum */
71 static ldrealtype VAllReduce_long_NrnParallelLD(ldrealtype d, int op, MPI_Comm comm);
72 #else
73 #define VAllReduce_long_NrnParallelLD VAllReduce_NrnParallelLD
74 #endif
75 /* z=x */
76 static void VCopy_NrnParallelLD(N_Vector x, N_Vector z);
77 /* z=x+y */
78 static void VSum_NrnParallelLD(N_Vector x, N_Vector y, N_Vector z);
79 /* z=x-y */
80 static void VDiff_NrnParallelLD(N_Vector x, N_Vector y, N_Vector z);
81 /* z=-x */
82 static void VNeg_NrnParallelLD(N_Vector x, N_Vector z);
83 /* z=c(x+y) */
84 static void VScaleSum_NrnParallelLD(realtype c, N_Vector x, N_Vector y, N_Vector z);
85 /* z=c(x-y) */
86 static void VScaleDiff_NrnParallelLD(realtype c, N_Vector x, N_Vector y, N_Vector z);
87 /* z=ax+y */
88 static void VLin1_NrnParallelLD(realtype a, N_Vector x, N_Vector y, N_Vector z);
89 /* z=ax-y */
90 static void VLin2_NrnParallelLD(realtype a, N_Vector x, N_Vector y, N_Vector z);
91 /* y <- ax+y */
92 static void Vaxpy_NrnParallelLD(realtype a, N_Vector x, N_Vector y);
93 /* x <- ax */
94 static void VScaleBy_NrnParallelLD(realtype a, N_Vector x);
95 
96 /*
97  * -----------------------------------------------------------------
98  * exported functions
99  * -----------------------------------------------------------------
100  */
101 
102 /* ----------------------------------------------------------------
103  * Function to create a new parallel vector with empty data array
104  */
105 
106 N_Vector N_VNewEmpty_NrnParallelLD(MPI_Comm comm, long int local_length, long int global_length) {
107  N_Vector v;
108  N_Vector_Ops ops;
110  long int n, Nsum;
111 
112  /* Compute global length as sum of local lengths */
113  n = local_length;
114 #if NRNMPI_DYNAMICLOAD
115  nrnmpi_long_allreduce_vec(&n, &Nsum, 1, 1);
116 #else
117  comm = nrnmpi_comm;
118  MPI_Allreduce(&n, &Nsum, 1, PVEC_INTEGER_MPI_TYPE, MPI_SUM, comm);
119 #endif
120  if (Nsum != global_length) {
121  printf(BAD_N);
122  return (NULL);
123  }
124 
125  /* Create vector */
126  v = (N_Vector) malloc(sizeof *v);
127  if (v == NULL)
128  return (NULL);
129 
130  /* Create vector operation structure */
131  ops = (N_Vector_Ops) malloc(sizeof(struct _generic_N_Vector_Ops));
132  if (ops == NULL) {
133  free(v);
134  return (NULL);
135  }
136 
137  ops->nvclone = N_VClone_NrnParallelLD;
138  ops->nvdestroy = N_VDestroy_NrnParallelLD;
139  ops->nvspace = N_VSpace_NrnParallelLD;
140  ops->nvgetarraypointer = N_VGetArrayPointer_NrnParallelLD;
141  ops->nvsetarraypointer = N_VSetArrayPointer_NrnParallelLD;
142  ops->nvlinearsum = N_VLinearSum_NrnParallelLD;
143  ops->nvconst = N_VConst_NrnParallelLD;
144  ops->nvprod = N_VProd_NrnParallelLD;
145  ops->nvdiv = N_VDiv_NrnParallelLD;
146  ops->nvscale = N_VScale_NrnParallelLD;
147  ops->nvabs = N_VAbs_NrnParallelLD;
148  ops->nvinv = N_VInv_NrnParallelLD;
149  ops->nvaddconst = N_VAddConst_NrnParallelLD;
150  ops->nvdotprod = N_VDotProd_NrnParallelLD;
151  ops->nvmaxnorm = N_VMaxNorm_NrnParallelLD;
152  ops->nvwrmsnormmask = N_VWrmsNormMask_NrnParallelLD;
153  ops->nvwrmsnorm = N_VWrmsNorm_NrnParallelLD;
154  ops->nvmin = N_VMin_NrnParallelLD;
155  ops->nvwl2norm = N_VWL2Norm_NrnParallelLD;
156  ops->nvl1norm = N_VL1Norm_NrnParallelLD;
157  ops->nvcompare = N_VCompare_NrnParallelLD;
158  ops->nvinvtest = N_VInvTest_NrnParallelLD;
159  ops->nvconstrmask = N_VConstrMask_NrnParallelLD;
160  ops->nvminquotient = N_VMinQuotient_NrnParallelLD;
161 
162  /* Create content */
163  content = (N_VectorContent_NrnParallelLD) malloc(sizeof(struct _N_VectorContent_NrnParallelLD));
164  if (content == NULL) {
165  free(ops);
166  free(v);
167  return (NULL);
168  }
169 
170  /* Attach lengths and communicator */
171  content->local_length = local_length;
172  content->global_length = global_length;
173  content->comm = comm;
174  content->own_data = FALSE;
175  content->data = NULL;
176 
177  /* Attach content and ops */
178  v->content = content;
179  v->ops = ops;
180 
181  return (v);
182 }
183 
184 /* ----------------------------------------------------------------
185  * Function to create a new parallel vector
186  */
187 
188 extern "C" N_Vector N_VNew_NrnParallelLD(MPI_Comm comm,
189  long int local_length,
190  long int global_length) {
191  N_Vector v;
192  realtype* data;
193 
194  v = N_VNewEmpty_NrnParallelLD(comm, local_length, global_length);
195  if (v == NULL)
196  return (NULL);
197 
198  /* Create data */
199  if (local_length > 0) {
200  /* Allocate memory */
201  data = (realtype*) malloc(local_length * sizeof(realtype));
202  if (data == NULL) {
204  return (NULL);
205  }
206 
207  /* Attach data */
209  NV_DATA_P_LD(v) = data;
210  }
211 
212  return (v);
213 }
214 
215 /* ----------------------------------------------------------------------------
216  * Function to clone from a template a new vector with empty (NULL) data array
217  */
218 
219 N_Vector N_VCloneEmpty_NrnParallelLD(N_Vector w) {
220  N_Vector v;
221  N_Vector_Ops ops;
223 
224  if (w == NULL)
225  return (NULL);
226 
227  /* Create vector */
228  v = (N_Vector) malloc(sizeof *v);
229  if (v == NULL)
230  return (NULL);
231 
232  /* Create vector operation structure */
233  ops = (N_Vector_Ops) malloc(sizeof(struct _generic_N_Vector_Ops));
234  if (ops == NULL) {
235  free(v);
236  return (NULL);
237  }
238 
239  ops->nvclone = w->ops->nvclone;
240  ops->nvdestroy = w->ops->nvdestroy;
241  ops->nvspace = w->ops->nvspace;
242  ops->nvgetarraypointer = w->ops->nvgetarraypointer;
243  ops->nvsetarraypointer = w->ops->nvsetarraypointer;
244  ops->nvlinearsum = w->ops->nvlinearsum;
245  ops->nvconst = w->ops->nvconst;
246  ops->nvprod = w->ops->nvprod;
247  ops->nvdiv = w->ops->nvdiv;
248  ops->nvscale = w->ops->nvscale;
249  ops->nvabs = w->ops->nvabs;
250  ops->nvinv = w->ops->nvinv;
251  ops->nvaddconst = w->ops->nvaddconst;
252  ops->nvdotprod = w->ops->nvdotprod;
253  ops->nvmaxnorm = w->ops->nvmaxnorm;
254  ops->nvwrmsnormmask = w->ops->nvwrmsnormmask;
255  ops->nvwrmsnorm = w->ops->nvwrmsnorm;
256  ops->nvmin = w->ops->nvmin;
257  ops->nvwl2norm = w->ops->nvwl2norm;
258  ops->nvl1norm = w->ops->nvl1norm;
259  ops->nvcompare = w->ops->nvcompare;
260  ops->nvinvtest = w->ops->nvinvtest;
261  ops->nvconstrmask = w->ops->nvconstrmask;
262  ops->nvminquotient = w->ops->nvminquotient;
263 
264  /* Create content */
265  content = (N_VectorContent_NrnParallelLD) malloc(sizeof(struct _N_VectorContent_NrnParallelLD));
266  if (content == NULL) {
267  free(ops);
268  free(v);
269  return (NULL);
270  }
271 
272  /* Attach lengths and communicator */
273  content->local_length = NV_LOCLENGTH_P_LD(w);
274  content->global_length = NV_GLOBLENGTH_P_LD(w);
275  content->comm = NV_COMM_P_LD(w);
276  content->own_data = FALSE;
277  content->data = NULL;
278 
279  /* Attach content and ops */
280  v->content = content;
281  v->ops = ops;
282 
283  return (v);
284 }
285 
286 /* ----------------------------------------------------------------
287  * Function to create a parallel N_Vector with user data component
288  */
289 
290 N_Vector N_VMake_NrnParallelLD(MPI_Comm comm,
291  long int local_length,
292  long int global_length,
293  realtype* v_data) {
294  N_Vector v;
295 
296  v = N_VNewEmpty_NrnParallelLD(comm, local_length, global_length);
297  if (v == NULL)
298  return (NULL);
299 
300  if (local_length > 0) {
301  /* Attach data */
303  NV_DATA_P_LD(v) = v_data;
304  }
305 
306  return (v);
307 }
308 
309 /* ----------------------------------------------------------------
310  * Function to create an array of new parallel vectors.
311  */
312 
314  MPI_Comm comm,
315  long int local_length,
316  long int global_length) {
317  N_Vector* vs;
318  int j;
319 
320  if (count <= 0)
321  return (NULL);
322 
323  vs = (N_Vector*) malloc(count * sizeof(N_Vector));
324  if (vs == NULL)
325  return (NULL);
326 
327  for (j = 0; j < count; j++) {
328  vs[j] = N_VNew_NrnParallelLD(comm, local_length, global_length);
329  if (vs[j] == NULL) {
331  return (NULL);
332  }
333  }
334 
335  return (vs);
336 }
337 
338 /* ----------------------------------------------------------------
339  * Function to create an array of new parallel vectors with empty
340  * (NULL) data array.
341  */
342 
344  MPI_Comm comm,
345  long int local_length,
346  long int global_length) {
347  N_Vector* vs;
348  int j;
349 
350  if (count <= 0)
351  return (NULL);
352 
353  vs = (N_Vector*) malloc(count * sizeof(N_Vector));
354  if (vs == NULL)
355  return (NULL);
356 
357  for (j = 0; j < count; j++) {
358  vs[j] = N_VNewEmpty_NrnParallelLD(comm, local_length, global_length);
359  if (vs[j] == NULL) {
361  return (NULL);
362  }
363  }
364 
365  return (vs);
366 }
367 
368 /* ----------------------------------------------------------------
369  * Function to free an array created with N_VNewVectorArray_NrnParallelLD
370  */
371 
372 void N_VDestroyVectorArray_NrnParallelLD(N_Vector* vs, int count) {
373  int j;
374 
375  for (j = 0; j < count; j++)
377 
378  free(vs);
379 }
380 
381 /* ----------------------------------------------------------------
382  * Function to print a parallel vector
383  */
384 
385 void N_VPrint_NrnParallelLD(N_Vector x) {
386  long int i, N;
387  realtype* xd;
388 
389  N = NV_LOCLENGTH_P_LD(x);
390  xd = NV_DATA_P_LD(x);
391 
392  for (i = 0; i < N; i++) {
393 #if defined(SUNDIALS_EXTENDED_PRECISION)
394  printf("%Lg\n", *xd++);
395 #elif defined(SUNDIALS_DOUBLE_PRECISION)
396  printf("%lg\n", *xd++);
397 #else
398  printf("%g\n", *xd++);
399 #endif
400  }
401  printf("\n");
402 }
403 
404 /*
405  * -----------------------------------------------------------------
406  * implementation of vector operations
407  * -----------------------------------------------------------------
408  */
409 
410 N_Vector N_VClone_NrnParallelLD(N_Vector w) {
411  N_Vector v;
412  realtype* data;
413  long int local_length;
414 
416  if (v == NULL)
417  return (NULL);
418 
419  local_length = NV_LOCLENGTH_P_LD(w);
420 
421  /* Create data */
422  if (local_length > 0) {
423  /* Allocate memory */
424  data = (realtype*) malloc(local_length * sizeof(realtype));
425  if (data == NULL) {
427  return (NULL);
428  }
429 
430  /* Attach data */
432  NV_DATA_P_LD(v) = data;
433  }
434 
435  return (v);
436 }
437 
438 void N_VDestroy_NrnParallelLD(N_Vector v) {
439  if ((NV_OWN_DATA_P_LD(v) == TRUE) && (NV_DATA_P_LD(v) != NULL))
440  free(NV_DATA_P_LD(v));
441  free(v->content);
442  free(v->ops);
443  free(v);
444 }
445 
446 void N_VSpace_NrnParallelLD(N_Vector v, long int* lrw, long int* liw) {
447  MPI_Comm comm;
448  int npes;
449 
450  comm = NV_COMM_P_LD(v);
451 #if NRNMPI_DYNAMICLOAD
452  npes = nrnmpi_numprocs;
453 #else
454  MPI_Comm_size(comm, &npes);
455 #endif
456 
457  *lrw = NV_GLOBLENGTH_P_LD(v);
458  *liw = 2 * npes;
459 }
460 
461 realtype* N_VGetArrayPointer_NrnParallelLD(N_Vector v) {
462  realtype* v_data;
463 
464  v_data = NV_DATA_P_LD(v);
465 
466  return (v_data);
467 }
468 
469 void N_VSetArrayPointer_NrnParallelLD(realtype* v_data, N_Vector v) {
470  if (NV_LOCLENGTH_P_LD(v) > 0)
471  NV_DATA_P_LD(v) = v_data;
472 }
473 
474 void N_VLinearSum_NrnParallelLD(realtype a, N_Vector x, realtype b, N_Vector y, N_Vector z) {
475  long int i, N;
476  realtype c, *xd, *yd, *zd;
477  N_Vector v1, v2;
478  booleantype test;
479 
480  if ((b == ONE) && (z == y)) { /* BLAS usage: axpy y <- ax+y */
481  Vaxpy_NrnParallelLD(a, x, y);
482  return;
483  }
484 
485  if ((a == ONE) && (z == x)) { /* BLAS usage: axpy x <- by+x */
486  Vaxpy_NrnParallelLD(b, y, x);
487  return;
488  }
489 
490  /* Case: a == b == 1.0 */
491 
492  if ((a == ONE) && (b == ONE)) {
493  VSum_NrnParallelLD(x, y, z);
494  return;
495  }
496 
497  /* Cases: (1) a == 1.0, b = -1.0, (2) a == -1.0, b == 1.0 */
498 
499  if ((test = ((a == ONE) && (b == -ONE))) || ((a == -ONE) && (b == ONE))) {
500  v1 = test ? y : x;
501  v2 = test ? x : y;
502  VDiff_NrnParallelLD(v2, v1, z);
503  return;
504  }
505 
506  /* Cases: (1) a == 1.0, b == other or 0.0, (2) a == other or 0.0, b == 1.0 */
507  /* if a or b is 0.0, then user should have called N_VScale */
508 
509  if ((test = (a == ONE)) || (b == ONE)) {
510  c = test ? b : a;
511  v1 = test ? y : x;
512  v2 = test ? x : y;
513  VLin1_NrnParallelLD(c, v1, v2, z);
514  return;
515  }
516 
517  /* Cases: (1) a == -1.0, b != 1.0, (2) a != 1.0, b == -1.0 */
518 
519  if ((test = (a == -ONE)) || (b == -ONE)) {
520  c = test ? b : a;
521  v1 = test ? y : x;
522  v2 = test ? x : y;
523  VLin2_NrnParallelLD(c, v1, v2, z);
524  return;
525  }
526 
527  /* Case: a == b */
528  /* catches case both a and b are 0.0 - user should have called N_VConst */
529 
530  if (a == b) {
531  VScaleSum_NrnParallelLD(a, x, y, z);
532  return;
533  }
534 
535  /* Case: a == -b */
536 
537  if (a == -b) {
538  VScaleDiff_NrnParallelLD(a, x, y, z);
539  return;
540  }
541 
542  /* Do all cases not handled above:
543  (1) a == other, b == 0.0 - user should have called N_VScale
544  (2) a == 0.0, b == other - user should have called N_VScale
545  (3) a,b == other, a !=b, a != -b */
546 
547  N = NV_LOCLENGTH_P_LD(x);
548  xd = NV_DATA_P_LD(x);
549  yd = NV_DATA_P_LD(y);
550  zd = NV_DATA_P_LD(z);
551 
552  for (i = 0; i < N; i++)
553  *zd++ = a * (*xd++) + b * (*yd++);
554 }
555 
556 void N_VConst_NrnParallelLD(realtype c, N_Vector z) {
557  long int i, N;
558  realtype* zd;
559 
560  N = NV_LOCLENGTH_P_LD(z);
561  zd = NV_DATA_P_LD(z);
562 
563  for (i = 0; i < N; i++)
564  *zd++ = c;
565 }
566 
567 void N_VProd_NrnParallelLD(N_Vector x, N_Vector y, N_Vector z) {
568  long int i, N;
569  realtype *xd, *yd, *zd;
570 
571  N = NV_LOCLENGTH_P_LD(x);
572  xd = NV_DATA_P_LD(x);
573  yd = NV_DATA_P_LD(y);
574  zd = NV_DATA_P_LD(z);
575 
576  for (i = 0; i < N; i++)
577  *zd++ = (*xd++) * (*yd++);
578 }
579 
580 void N_VDiv_NrnParallelLD(N_Vector x, N_Vector y, N_Vector z) {
581  long int i, N;
582  realtype *xd, *yd, *zd;
583 
584  N = NV_LOCLENGTH_P_LD(x);
585  xd = NV_DATA_P_LD(x);
586  yd = NV_DATA_P_LD(y);
587  zd = NV_DATA_P_LD(z);
588 
589  for (i = 0; i < N; i++)
590  *zd++ = (*xd++) / (*yd++);
591 }
592 
593 void N_VScale_NrnParallelLD(realtype c, N_Vector x, N_Vector z) {
594  long int i, N;
595  realtype *xd, *zd;
596 
597  if (z == x) { /* BLAS usage: scale x <- cx */
599  return;
600  }
601 
602  if (c == ONE) {
603  VCopy_NrnParallelLD(x, z);
604  } else if (c == -ONE) {
605  VNeg_NrnParallelLD(x, z);
606  } else {
607  N = NV_LOCLENGTH_P_LD(x);
608  xd = NV_DATA_P_LD(x);
609  zd = NV_DATA_P_LD(z);
610  for (i = 0; i < N; i++)
611  *zd++ = c * (*xd++);
612  }
613 }
614 
615 void N_VAbs_NrnParallelLD(N_Vector x, N_Vector z) {
616  long int i, N;
617  realtype *xd, *zd;
618 
619  N = NV_LOCLENGTH_P_LD(x);
620  xd = NV_DATA_P_LD(x);
621  zd = NV_DATA_P_LD(z);
622 
623  for (i = 0; i < N; i++, xd++, zd++)
624  *zd = ABS(*xd);
625 }
626 
627 void N_VInv_NrnParallelLD(N_Vector x, N_Vector z) {
628  long int i, N;
629  realtype *xd, *zd;
630 
631  N = NV_LOCLENGTH_P_LD(x);
632  xd = NV_DATA_P_LD(x);
633  zd = NV_DATA_P_LD(z);
634 
635  for (i = 0; i < N; i++)
636  *zd++ = ONE / (*xd++);
637 }
638 
639 void N_VAddConst_NrnParallelLD(N_Vector x, realtype b, N_Vector z) {
640  long int i, N;
641  realtype *xd, *zd;
642 
643  N = NV_LOCLENGTH_P_LD(x);
644  xd = NV_DATA_P_LD(x);
645  zd = NV_DATA_P_LD(z);
646 
647  for (i = 0; i < N; i++)
648  *zd++ = (*xd++) + b;
649 }
650 
651 realtype N_VDotProd_NrnParallelLD(N_Vector x, N_Vector y) {
652  long int i, N;
653  realtype sum = ZERO, *xd, *yd, gsum;
654  MPI_Comm comm;
655 
656  N = NV_LOCLENGTH_P_LD(x);
657  xd = NV_DATA_P_LD(x);
658  yd = NV_DATA_P_LD(y);
659  comm = NV_COMM_P_LD(x);
660 
661  for (i = 0; i < N; i++)
662  sum += xd[i] * yd[i];
663 
664  gsum = VAllReduce_long_NrnParallelLD(sum, 1, comm);
665  return (gsum);
666 }
667 
668 realtype N_VMaxNorm_NrnParallelLD(N_Vector x) {
669  long int i, N;
670  realtype max, *xd, gmax;
671  MPI_Comm comm;
672 
673  N = NV_LOCLENGTH_P_LD(x);
674  xd = NV_DATA_P_LD(x);
675  comm = NV_COMM_P_LD(x);
676 
677  max = ZERO;
678 
679  for (i = 0; i < N; i++, xd++) {
680  if (ABS(*xd) > max)
681  max = ABS(*xd);
682  }
683 
684  gmax = VAllReduce_NrnParallelLD(max, 2, comm);
685  return (gmax);
686 }
687 
688 realtype N_VWrmsNorm_NrnParallelLD(N_Vector x, N_Vector w) {
689  long int i, N, N_global;
690  realtype prodi, *xd, *wd;
691  ldrealtype sum = ZERO, gsum;
692  MPI_Comm comm;
693 
694  N = NV_LOCLENGTH_P_LD(x);
695  N_global = NV_GLOBLENGTH_P_LD(x);
696  xd = NV_DATA_P_LD(x);
697  wd = NV_DATA_P_LD(w);
698  comm = NV_COMM_P_LD(x);
699 
700  for (i = 0; i < N; i++) {
701  prodi = (*xd++) * (*wd++);
702  sum += prodi * prodi;
703  }
704 
705  gsum = VAllReduce_long_NrnParallelLD(sum, 1, comm);
706  return (RSqrt((realtype) gsum / N_global));
707 }
708 
709 realtype N_VWrmsNormMask_NrnParallelLD(N_Vector x, N_Vector w, N_Vector id) {
710  long int i, N, N_global;
711  realtype prodi, *xd, *wd, *idd;
712  ldrealtype sum = ZERO, gsum;
713  MPI_Comm comm;
714 
715  N = NV_LOCLENGTH_P_LD(x);
716  N_global = NV_GLOBLENGTH_P_LD(x);
717  xd = NV_DATA_P_LD(x);
718  wd = NV_DATA_P_LD(w);
719  idd = NV_DATA_P_LD(id);
720  comm = NV_COMM_P_LD(x);
721 
722  for (i = 0; i < N; i++) {
723  if (idd[i] > ZERO) {
724  prodi = xd[i] * wd[i];
725  sum += prodi * prodi;
726  }
727  }
728 
729  gsum = VAllReduce_long_NrnParallelLD(sum, 1, comm);
730  return (RSqrt((realtype) gsum / N_global));
731 }
732 
733 realtype N_VMin_NrnParallelLD(N_Vector x) {
734  long int i, N;
735  realtype min, *xd, gmin;
736  MPI_Comm comm;
737 
738  N = NV_LOCLENGTH_P_LD(x);
739  comm = NV_COMM_P_LD(x);
740 
741  min = BIG_REAL;
742 
743  if (N > 0) {
744  xd = NV_DATA_P_LD(x);
745 
746  min = xd[0];
747 
748  xd++;
749  for (i = 1; i < N; i++, xd++) {
750  if ((*xd) < min)
751  min = *xd;
752  }
753  }
754 
755  gmin = VAllReduce_NrnParallelLD(min, 3, comm);
756  return (gmin);
757 }
758 
759 realtype N_VWL2Norm_NrnParallelLD(N_Vector x, N_Vector w) {
760  long int i, N;
761  realtype prodi, *xd, *wd;
762  ldrealtype sum = ZERO, gsum;
763  MPI_Comm comm;
764 
765  N = NV_LOCLENGTH_P_LD(x);
766  xd = NV_DATA_P_LD(x);
767  wd = NV_DATA_P_LD(w);
768  comm = NV_COMM_P_LD(x);
769 
770  for (i = 0; i < N; i++) {
771  prodi = (*xd++) * (*wd++);
772  sum += prodi * prodi;
773  }
774 
775  gsum = VAllReduce_long_NrnParallelLD(sum, 1, comm);
776  return (RSqrt((realtype) gsum));
777 }
778 
779 realtype N_VL1Norm_NrnParallelLD(N_Vector x) {
780  long int i, N;
781  realtype* xd;
782  ldrealtype sum = ZERO, gsum;
783  MPI_Comm comm;
784 
785  N = NV_LOCLENGTH_P_LD(x);
786  xd = NV_DATA_P_LD(x);
787  comm = NV_COMM_P_LD(x);
788 
789  for (i = 0; i < N; i++, xd++)
790  sum += ABS(*xd);
791 
792  gsum = VAllReduce_long_NrnParallelLD((realtype) sum, 1, comm);
793  return (gsum);
794 }
795 
796 void N_VCompare_NrnParallelLD(realtype c, N_Vector x, N_Vector z) {
797  long int i, N;
798  realtype *xd, *zd;
799 
800  N = NV_LOCLENGTH_P_LD(x);
801  xd = NV_DATA_P_LD(x);
802  zd = NV_DATA_P_LD(z);
803 
804  for (i = 0; i < N; i++, xd++, zd++) {
805  *zd = (ABS(*xd) >= c) ? ONE : ZERO;
806  }
807 }
808 
809 booleantype N_VInvTest_NrnParallelLD(N_Vector x, N_Vector z) {
810  long int i, N;
811  realtype *xd, *zd, val, gval;
812  MPI_Comm comm;
813 
814  N = NV_LOCLENGTH_P_LD(x);
815  xd = NV_DATA_P_LD(x);
816  zd = NV_DATA_P_LD(z);
817  comm = NV_COMM_P_LD(x);
818 
819  val = ONE;
820  for (i = 0; i < N; i++) {
821  if (*xd == ZERO)
822  val = ZERO;
823  else
824  *zd++ = ONE / (*xd++);
825  }
826 
827  gval = VAllReduce_NrnParallelLD(val, 3, comm);
828  if (gval == ZERO)
829  return (FALSE);
830  else
831  return (TRUE);
832 }
833 
834 booleantype N_VConstrMask_NrnParallelLD(N_Vector c, N_Vector x, N_Vector m) {
835  long int i, N;
836  booleantype test;
837  realtype *cd, *xd, *md;
838  MPI_Comm comm;
839 
840  N = NV_LOCLENGTH_P_LD(x);
841  xd = NV_DATA_P_LD(x);
842  cd = NV_DATA_P_LD(c);
843  md = NV_DATA_P_LD(m);
844  comm = NV_COMM_P_LD(x);
845 
846  test = TRUE;
847 
848  for (i = 0; i < N; i++, cd++, xd++, md++) {
849  *md = ZERO;
850  if (*cd == ZERO)
851  continue;
852  if (*cd > ONEPT5 || (*cd) < -ONEPT5) {
853  if ((*xd) * (*cd) <= ZERO) {
854  test = FALSE;
855  *md = ONE;
856  }
857  continue;
858  }
859  if ((*cd) > HALF || (*cd) < -HALF) {
860  if ((*xd) * (*cd) < ZERO) {
861  test = FALSE;
862  *md = ONE;
863  }
864  }
865  }
866 
867  return ((booleantype) VAllReduce_long_NrnParallelLD((realtype) test, 3, comm));
868 }
869 
870 realtype N_VMinQuotient_NrnParallelLD(N_Vector num, N_Vector denom) {
871  booleantype notEvenOnce;
872  long int i, N;
873  realtype *nd, *dd, min = 0.0;
874  MPI_Comm comm;
875 
876  N = NV_LOCLENGTH_P_LD(num);
877  nd = NV_DATA_P_LD(num);
878  dd = NV_DATA_P_LD(denom);
879  comm = NV_COMM_P_LD(num);
880 
881  notEvenOnce = TRUE;
882 
883  for (i = 0; i < N; i++, nd++, dd++) {
884  if (*dd == ZERO)
885  continue;
886  else {
887  if (notEvenOnce) {
888  min = *nd / *dd;
889  notEvenOnce = FALSE;
890  } else
891  min = MIN(min, (*nd) / (*dd));
892  }
893  }
894 
895  if (notEvenOnce || (N == 0))
896  min = BIG_REAL;
897 
898  return (VAllReduce_NrnParallelLD(min, 3, comm));
899 }
900 
901 /*
902  * -----------------------------------------------------------------
903  * private functions
904  * -----------------------------------------------------------------
905  */
906 
907 static realtype VAllReduce_NrnParallelLD(realtype d, int op, MPI_Comm comm) {
908  /*
909  * This function does a global reduction. The operation is
910  * sum if op = 1,
911  * max if op = 2,
912  * min if op = 3.
913  * The operation is over all processors in the communicator
914  */
915 
916  realtype out;
917 
918 #if NRNMPI_DYNAMICLOAD
919  nrnmpi_dbl_allreduce_vec(&d, &out, 1, op);
920 #else
921  switch (op) {
922  case 1:
923  MPI_Allreduce(&d, &out, 1, PVEC_REAL_MPI_TYPE, MPI_SUM, comm);
924  break;
925 
926  case 2:
927  MPI_Allreduce(&d, &out, 1, PVEC_REAL_MPI_TYPE, MPI_MAX, comm);
928  break;
929 
930  case 3:
931  MPI_Allreduce(&d, &out, 1, PVEC_REAL_MPI_TYPE, MPI_MIN, comm);
932  break;
933 
934  default:
935  break;
936  }
937 #endif
938 
939  return (out);
940 }
941 
942 #if USELONGDOUBLE
943 static ldrealtype VAllReduce_long_NrnParallelLD(ldrealtype d, int op, MPI_Comm comm) {
944  /*
945  * This function does a global reduction. The operation is
946  * sum if op = 1,
947  * max if op = 2,
948  * min if op = 3.
949  * The operation is over all processors in the communicator
950  */
951 
952  ldrealtype out;
953 
954  assert(op == 1);
955 #if NRNMPI_DYNAMICLOAD
956  nrnmpi_longdbl_allreduce_vec(&d, &out, 1, op);
957 #else
958  MPI_Allreduce(&d, &out, 1, MPI_LONG_DOUBLE, MPI_SUM, comm);
959 #endif
960 
961  return (out);
962 }
963 #endif
964 
965 static void VCopy_NrnParallelLD(N_Vector x, N_Vector z) {
966  long int i, N;
967  realtype *xd, *zd;
968 
969  N = NV_LOCLENGTH_P_LD(x);
970  xd = NV_DATA_P_LD(x);
971  zd = NV_DATA_P_LD(z);
972 
973  for (i = 0; i < N; i++)
974  *zd++ = *xd++;
975 }
976 
977 static void VSum_NrnParallelLD(N_Vector x, N_Vector y, N_Vector z) {
978  long int i, N;
979  realtype *xd, *yd, *zd;
980 
981  N = NV_LOCLENGTH_P_LD(x);
982  xd = NV_DATA_P_LD(x);
983  yd = NV_DATA_P_LD(y);
984  zd = NV_DATA_P_LD(z);
985 
986  for (i = 0; i < N; i++)
987  *zd++ = (*xd++) + (*yd++);
988 }
989 
990 static void VDiff_NrnParallelLD(N_Vector x, N_Vector y, N_Vector z) {
991  long int i, N;
992  realtype *xd, *yd, *zd;
993 
994  N = NV_LOCLENGTH_P_LD(x);
995  xd = NV_DATA_P_LD(x);
996  yd = NV_DATA_P_LD(y);
997  zd = NV_DATA_P_LD(z);
998 
999  for (i = 0; i < N; i++)
1000  *zd++ = (*xd++) - (*yd++);
1001 }
1002 
1003 static void VNeg_NrnParallelLD(N_Vector x, N_Vector z) {
1004  long int i, N;
1005  realtype *xd, *zd;
1006 
1007  N = NV_LOCLENGTH_P_LD(x);
1008  xd = NV_DATA_P_LD(x);
1009  zd = NV_DATA_P_LD(z);
1010 
1011  for (i = 0; i < N; i++)
1012  *zd++ = -(*xd++);
1013 }
1014 
1015 static void VScaleSum_NrnParallelLD(realtype c, N_Vector x, N_Vector y, N_Vector z) {
1016  long int i, N;
1017  realtype *xd, *yd, *zd;
1018 
1019  N = NV_LOCLENGTH_P_LD(x);
1020  xd = NV_DATA_P_LD(x);
1021  yd = NV_DATA_P_LD(y);
1022  zd = NV_DATA_P_LD(z);
1023 
1024  for (i = 0; i < N; i++)
1025  *zd++ = c * ((*xd++) + (*yd++));
1026 }
1027 
1028 static void VScaleDiff_NrnParallelLD(realtype c, N_Vector x, N_Vector y, N_Vector z) {
1029  long int i, N;
1030  realtype *xd, *yd, *zd;
1031 
1032  N = NV_LOCLENGTH_P_LD(x);
1033  xd = NV_DATA_P_LD(x);
1034  yd = NV_DATA_P_LD(y);
1035  zd = NV_DATA_P_LD(z);
1036 
1037  for (i = 0; i < N; i++)
1038  *zd++ = c * ((*xd++) - (*yd++));
1039 }
1040 
1041 static void VLin1_NrnParallelLD(realtype a, N_Vector x, N_Vector y, N_Vector z) {
1042  long int i, N;
1043  realtype *xd, *yd, *zd;
1044 
1045  N = NV_LOCLENGTH_P_LD(x);
1046  xd = NV_DATA_P_LD(x);
1047  yd = NV_DATA_P_LD(y);
1048  zd = NV_DATA_P_LD(z);
1049 
1050  for (i = 0; i < N; i++)
1051  *zd++ = a * (*xd++) + (*yd++);
1052 }
1053 
1054 static void VLin2_NrnParallelLD(realtype a, N_Vector x, N_Vector y, N_Vector z) {
1055  long int i, N;
1056  realtype *xd, *yd, *zd;
1057 
1058  N = NV_LOCLENGTH_P_LD(x);
1059  xd = NV_DATA_P_LD(x);
1060  yd = NV_DATA_P_LD(y);
1061  zd = NV_DATA_P_LD(z);
1062 
1063  for (i = 0; i < N; i++)
1064  *zd++ = a * (*xd++) - (*yd++);
1065 }
1066 
1067 static void Vaxpy_NrnParallelLD(realtype a, N_Vector x, N_Vector y) {
1068  long int i, N;
1069  realtype *xd, *yd;
1070 
1071  N = NV_LOCLENGTH_P_LD(x);
1072  xd = NV_DATA_P_LD(x);
1073  yd = NV_DATA_P_LD(y);
1074 
1075  if (a == ONE) {
1076  for (i = 0; i < N; i++)
1077  *yd++ += (*xd++);
1078  return;
1079  }
1080 
1081  if (a == -ONE) {
1082  for (i = 0; i < N; i++)
1083  *yd++ -= (*xd++);
1084  return;
1085  }
1086 
1087  for (i = 0; i < N; i++)
1088  *yd++ += a * (*xd++);
1089 }
1090 
1091 static void VScaleBy_NrnParallelLD(realtype a, N_Vector x) {
1092  long int i, N;
1093  realtype* xd;
1094 
1095  N = NV_LOCLENGTH_P_LD(x);
1096  xd = NV_DATA_P_LD(x);
1097 
1098  for (i = 0; i < N; i++)
1099  *xd++ *= a;
1100 }
short type
Definition: cabvars.h:9
#define c
#define TRUE
Definition: err.c:57
#define FALSE
Definition: err.c:56
#define ABS(x)
Definition: extras.c:201
#define MIN(a, b)
Definition: grids.h:36
#define assert(ex)
Definition: hocassrt.h:32
#define min(a, b)
Definition: matrix.h:157
#define max(a, b)
Definition: matrix.h:154
#define v
Definition: md1redef.h:4
#define i
Definition: md1redef.h:12
#define printf
Definition: mwprefix.h:26
int const size_t const size_t n
Definition: nrngsl.h:11
size_t j
int nrnmpi_numprocs
static void VDiff_NrnParallelLD(N_Vector x, N_Vector y, N_Vector z)
static realtype VAllReduce_NrnParallelLD(realtype d, int op, MPI_Comm comm)
realtype N_VWL2Norm_NrnParallelLD(N_Vector x, N_Vector w)
void N_VSetArrayPointer_NrnParallelLD(realtype *v_data, N_Vector v)
N_Vector N_VClone_NrnParallelLD(N_Vector w)
booleantype N_VInvTest_NrnParallelLD(N_Vector x, N_Vector z)
static void VSum_NrnParallelLD(N_Vector x, N_Vector y, N_Vector z)
realtype N_VL1Norm_NrnParallelLD(N_Vector x)
void N_VProd_NrnParallelLD(N_Vector x, N_Vector y, N_Vector z)
#define ONE
void N_VScale_NrnParallelLD(realtype c, N_Vector x, N_Vector z)
static void VCopy_NrnParallelLD(N_Vector x, N_Vector z)
N_Vector N_VCloneEmpty_NrnParallelLD(N_Vector w)
N_Vector N_VMake_NrnParallelLD(MPI_Comm comm, long int local_length, long int global_length, realtype *v_data)
realtype N_VMin_NrnParallelLD(N_Vector x)
#define BAD_N
#define HALF
realtype N_VWrmsNormMask_NrnParallelLD(N_Vector x, N_Vector w, N_Vector id)
void N_VDiv_NrnParallelLD(N_Vector x, N_Vector y, N_Vector z)
void N_VSpace_NrnParallelLD(N_Vector v, long int *lrw, long int *liw)
N_Vector N_VNewEmpty_NrnParallelLD(MPI_Comm comm, long int local_length, long int global_length)
void N_VAddConst_NrnParallelLD(N_Vector x, realtype b, N_Vector z)
static void VScaleDiff_NrnParallelLD(realtype c, N_Vector x, N_Vector y, N_Vector z)
booleantype N_VConstrMask_NrnParallelLD(N_Vector c, N_Vector x, N_Vector m)
void N_VInv_NrnParallelLD(N_Vector x, N_Vector z)
void N_VAbs_NrnParallelLD(N_Vector x, N_Vector z)
static void VScaleBy_NrnParallelLD(realtype a, N_Vector x)
static ldrealtype VAllReduce_long_NrnParallelLD(ldrealtype d, int op, MPI_Comm comm)
void N_VLinearSum_NrnParallelLD(realtype a, N_Vector x, realtype b, N_Vector y, N_Vector z)
void N_VConst_NrnParallelLD(realtype c, N_Vector z)
void N_VDestroyVectorArray_NrnParallelLD(N_Vector *vs, int count)
void N_VPrint_NrnParallelLD(N_Vector x)
static void VNeg_NrnParallelLD(N_Vector x, N_Vector z)
static void VLin1_NrnParallelLD(realtype a, N_Vector x, N_Vector y, N_Vector z)
void N_VDestroy_NrnParallelLD(N_Vector v)
MPI_Comm nrnmpi_comm
#define ldrealtype
realtype N_VMaxNorm_NrnParallelLD(N_Vector x)
N_Vector * N_VNewVectorArray_NrnParallelLD(int count, MPI_Comm comm, long int local_length, long int global_length)
realtype * N_VGetArrayPointer_NrnParallelLD(N_Vector v)
realtype N_VDotProd_NrnParallelLD(N_Vector x, N_Vector y)
realtype N_VMinQuotient_NrnParallelLD(N_Vector num, N_Vector denom)
#define ZERO
realtype N_VWrmsNorm_NrnParallelLD(N_Vector x, N_Vector w)
void N_VCompare_NrnParallelLD(realtype c, N_Vector x, N_Vector z)
N_Vector N_VNew_NrnParallelLD(MPI_Comm comm, long int local_length, long int global_length)
#define ONEPT5
static void VScaleSum_NrnParallelLD(realtype c, N_Vector x, N_Vector y, N_Vector z)
N_Vector * N_VNewVectorArrayEmpty_NrnParallelLD(int count, MPI_Comm comm, long int local_length, long int global_length)
static void VLin2_NrnParallelLD(realtype a, N_Vector x, N_Vector y, N_Vector z)
static void Vaxpy_NrnParallelLD(realtype a, N_Vector x, N_Vector y)
#define NV_OWN_DATA_P_LD(v)
#define PVEC_INTEGER_MPI_TYPE
#define NV_COMM_P_LD(v)
struct _N_VectorContent_NrnParallelLD * N_VectorContent_NrnParallelLD
#define NV_DATA_P_LD(v)
#define NV_GLOBLENGTH_P_LD(v)
#define NV_LOCLENGTH_P_LD(v)
#define data
Definition: rbtqueue.cpp:49
#define cnt
Definition: spt2queue.cpp:19
#define NULL
Definition: sptree.h:16