19 #define USELONGDOUBLE 1
21 #include <../../nrnconf.h>
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,
35 extern "C" void nrnmpi_long_allreduce_vec(
long* src,
long* dest,
int cnt,
int type);
40 #if NRNMPI_DYNAMICLOAD
44 #include "sundialsmath.h"
45 #include "sundialstypes.h"
47 #define ZERO RCONST(0.0)
48 #define HALF RCONST(0.5)
49 #define ONE RCONST(1.0)
50 #define ONEPT5 RCONST(1.5)
53 #define ldrealtype long double
55 #define ldrealtype realtype
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
73 #define VAllReduce_long_NrnParallelLD VAllReduce_NrnParallelLD
114 #if NRNMPI_DYNAMICLOAD
115 nrnmpi_long_allreduce_vec(&
n, &Nsum, 1, 1);
120 if (Nsum != global_length) {
126 v = (N_Vector) malloc(
sizeof *
v);
131 ops = (N_Vector_Ops) malloc(
sizeof(
struct _generic_N_Vector_Ops));
164 if (content ==
NULL) {
173 content->
comm = comm;
178 v->content = content;
189 long int local_length,
190 long int global_length) {
199 if (local_length > 0) {
201 data = (realtype*) malloc(local_length *
sizeof(realtype));
228 v = (N_Vector) malloc(
sizeof *
v);
233 ops = (N_Vector_Ops) malloc(
sizeof(
struct _generic_N_Vector_Ops));
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;
266 if (content ==
NULL) {
280 v->content = content;
291 long int local_length,
292 long int global_length,
300 if (local_length > 0) {
315 long int local_length,
316 long int global_length) {
323 vs = (N_Vector*) malloc(count *
sizeof(N_Vector));
327 for (
j = 0;
j < count;
j++) {
345 long int local_length,
346 long int global_length) {
353 vs = (N_Vector*) malloc(count *
sizeof(N_Vector));
357 for (
j = 0;
j < count;
j++) {
375 for (
j = 0;
j < count;
j++)
392 for (
i = 0;
i < N;
i++) {
393 #if defined(SUNDIALS_EXTENDED_PRECISION)
395 #elif defined(SUNDIALS_DOUBLE_PRECISION)
413 long int local_length;
422 if (local_length > 0) {
424 data = (realtype*) malloc(local_length *
sizeof(realtype));
451 #if NRNMPI_DYNAMICLOAD
454 MPI_Comm_size(comm, &npes);
476 realtype
c, *xd, *yd, *zd;
480 if ((b ==
ONE) && (z == y)) {
485 if ((a ==
ONE) && (z == x)) {
492 if ((a ==
ONE) && (b ==
ONE)) {
499 if ((test = ((a ==
ONE) && (b == -
ONE))) || ((a == -
ONE) && (b ==
ONE))) {
509 if ((test = (a ==
ONE)) || (b ==
ONE)) {
519 if ((test = (a == -
ONE)) || (b == -
ONE)) {
552 for (
i = 0;
i < N;
i++)
553 *zd++ = a * (*xd++) + b * (*yd++);
563 for (
i = 0;
i < N;
i++)
569 realtype *xd, *yd, *zd;
576 for (
i = 0;
i < N;
i++)
577 *zd++ = (*xd++) * (*yd++);
582 realtype *xd, *yd, *zd;
589 for (
i = 0;
i < N;
i++)
590 *zd++ = (*xd++) / (*yd++);
604 }
else if (
c == -
ONE) {
610 for (
i = 0;
i < N;
i++)
623 for (
i = 0;
i < N;
i++, xd++, zd++)
635 for (
i = 0;
i < N;
i++)
636 *zd++ =
ONE / (*xd++);
647 for (
i = 0;
i < N;
i++)
653 realtype sum =
ZERO, *xd, *yd, gsum;
661 for (
i = 0;
i < N;
i++)
662 sum += xd[
i] * yd[
i];
670 realtype
max, *xd, gmax;
679 for (
i = 0;
i < N;
i++, xd++) {
689 long int i, N, N_global;
690 realtype prodi, *xd, *wd;
700 for (
i = 0;
i < N;
i++) {
701 prodi = (*xd++) * (*wd++);
702 sum += prodi * prodi;
706 return (RSqrt((realtype) gsum / N_global));
710 long int i, N, N_global;
711 realtype prodi, *xd, *wd, *idd;
722 for (
i = 0;
i < N;
i++) {
724 prodi = xd[
i] * wd[
i];
725 sum += prodi * prodi;
730 return (RSqrt((realtype) gsum / N_global));
735 realtype
min, *xd, gmin;
749 for (
i = 1;
i < N;
i++, xd++) {
761 realtype prodi, *xd, *wd;
770 for (
i = 0;
i < N;
i++) {
771 prodi = (*xd++) * (*wd++);
772 sum += prodi * prodi;
776 return (RSqrt((realtype) gsum));
789 for (
i = 0;
i < N;
i++, xd++)
804 for (
i = 0;
i < N;
i++, xd++, zd++) {
811 realtype *xd, *zd, val, gval;
820 for (
i = 0;
i < N;
i++) {
824 *zd++ =
ONE / (*xd++);
837 realtype *cd, *xd, *md;
848 for (
i = 0;
i < N;
i++, cd++, xd++, md++) {
853 if ((*xd) * (*cd) <=
ZERO) {
859 if ((*cd) >
HALF || (*cd) < -
HALF) {
860 if ((*xd) * (*cd) <
ZERO) {
871 booleantype notEvenOnce;
873 realtype *nd, *dd,
min = 0.0;
883 for (
i = 0;
i < N;
i++, nd++, dd++) {
895 if (notEvenOnce || (N == 0))
918 #if NRNMPI_DYNAMICLOAD
919 nrnmpi_dbl_allreduce_vec(&d, &out, 1, op);
923 MPI_Allreduce(&d, &out, 1, PVEC_REAL_MPI_TYPE, MPI_SUM, comm);
927 MPI_Allreduce(&d, &out, 1, PVEC_REAL_MPI_TYPE, MPI_MAX, comm);
931 MPI_Allreduce(&d, &out, 1, PVEC_REAL_MPI_TYPE, MPI_MIN, comm);
955 #if NRNMPI_DYNAMICLOAD
956 nrnmpi_longdbl_allreduce_vec(&d, &out, 1, op);
958 MPI_Allreduce(&d, &out, 1, MPI_LONG_DOUBLE, MPI_SUM, comm);
973 for (
i = 0;
i < N;
i++)
979 realtype *xd, *yd, *zd;
986 for (
i = 0;
i < N;
i++)
987 *zd++ = (*xd++) + (*yd++);
992 realtype *xd, *yd, *zd;
999 for (
i = 0;
i < N;
i++)
1000 *zd++ = (*xd++) - (*yd++);
1011 for (
i = 0;
i < N;
i++)
1017 realtype *xd, *yd, *zd;
1024 for (
i = 0;
i < N;
i++)
1025 *zd++ =
c * ((*xd++) + (*yd++));
1030 realtype *xd, *yd, *zd;
1037 for (
i = 0;
i < N;
i++)
1038 *zd++ =
c * ((*xd++) - (*yd++));
1043 realtype *xd, *yd, *zd;
1050 for (
i = 0;
i < N;
i++)
1051 *zd++ = a * (*xd++) + (*yd++);
1056 realtype *xd, *yd, *zd;
1063 for (
i = 0;
i < N;
i++)
1064 *zd++ = a * (*xd++) - (*yd++);
1076 for (
i = 0;
i < N;
i++)
1082 for (
i = 0;
i < N;
i++)
1087 for (
i = 0;
i < N;
i++)
1088 *yd++ += a * (*xd++);
1098 for (
i = 0;
i < N;
i++)
int const size_t const size_t n
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)
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)
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)
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)
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)
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
struct _N_VectorContent_NrnParallelLD * N_VectorContent_NrnParallelLD
#define NV_GLOBLENGTH_P_LD(v)
#define NV_LOCLENGTH_P_LD(v)