jacobi.template
Loading...
Searching...
No Matches
vdop.cpp
Go to the documentation of this file.
1#include "vdop.h"
2#include <cassert> // assert()
3#include <cmath>
4#include <iostream>
5#include <vector>
6using namespace std;
7
8
9void vddiv(vector<double> & x, vector<double> const& y,
10 vector<double> const& z)
11{
12 assert( x.size()==y.size() && y.size()==z.size() );
13 size_t n = x.size();
14
15#pragma omp parallel for
16 for (size_t k = 0; k < n; ++k)
17 {
18 x[k] = y[k] / z[k];
19 }
20 return;
21}
22
23//******************************************************************************
24
25void vdaxpy(std::vector<double> & x, std::vector<double> const& y,
26 double alpha, std::vector<double> const& z )
27{
28 assert( x.size()==y.size() && y.size()==z.size() );
29 size_t n = x.size();
30
31#pragma omp parallel for
32 for (size_t k = 0; k < n; ++k)
33 {
34 x[k] = y[k] + alpha * z[k];
35 }
36 return;
37}
38//******************************************************************************
39
40double dscapr(std::vector<double> const& x, std::vector<double> const& y)
41{
42 assert( x.size()==y.size());
43 size_t n = x.size();
44
45 double s = 0.0;
46#pragma omp parallel for reduction(+:s)
47 for (size_t k = 0; k < n; ++k)
48 {
49 s += x[k] * y[k];
50 }
51
52 return s;
53}
54
55//******************************************************************************
56bool CompareVectors(std::vector<double> const& x, int const n, double const y[], double const eps)
57{
58 bool bn = (static_cast<int>(x.size())==n);
59 if (!bn)
60 {
61 cout << "######### Error: " << "number of elements" << endl;
62 }
63 //bool bv = equal(x.cbegin(),x.cend(),y);
64 bool bv = equal(x.cbegin(),x.cend(),y,
65 [eps](double a, double b) -> bool
66 { return std::abs(a-b)<eps*(1.0+0.5*(std::abs(a)+ std::abs(b))); }
67 );
68 if (!bv)
69 {
70 assert(static_cast<int>(x.size())==n);
71 cout << "######### Error: " << "values" << endl;
72 }
73 return bn && bv;
74}
75
76//******************************************************************************
77double par_scalar(vector<double> const &x, vector<double> const &y, MPI_Comm const& icomm)
78{
79 const double s = dscapr(x,y);
80 double sg;
81 MPI_Allreduce(&s,&sg,1,MPI_DOUBLE,MPI_SUM,icomm);
82
83 return(sg);
84}
85
86//******************************************************************************
87void ExchangeAll(vector<double> const &xin, vector<double> &yout, MPI_Comm const &icomm)
88{
89 int myrank, numprocs,ierr(-1);
90 MPI_Comm_rank(icomm, &myrank); // my MPI-rank
91 MPI_Comm_size(icomm, &numprocs);
92 int const N=xin.size();
93 int const sendcount = N/numprocs; // equal sized junks
94 assert(sendcount*numprocs==N); // really all junk sized?
95 assert(xin.size()==yout.size());
96
97 auto sendbuf = xin.data();
98 auto recvbuf = yout.data();
99 ierr = MPI_Alltoall(sendbuf, sendcount, MPI_DOUBLE,
100 recvbuf, sendcount, MPI_DOUBLE, icomm);
101 assert(0==ierr);
102
103 return;
104}
105
106//******************************************************************************
107void ExchangeAllInPlace(vector<double> &xin, MPI_Comm const &icomm)
108{
109 int myrank, numprocs,ierr(-1);
110 MPI_Comm_rank(icomm, &myrank); // my MPI-rank
111 MPI_Comm_size(icomm, &numprocs);
112 int const N=xin.size();
113 int const sendcount = N/numprocs; // equal sized junks
114 assert(sendcount*numprocs==N); // really all junk sized?
115
116 auto sendbuf = xin.data();
117 ierr = MPI_Alltoall(MPI_IN_PLACE, sendcount, MPI_DOUBLE,
118 sendbuf, sendcount, MPI_DOUBLE, icomm);
119 assert(0==ierr);
120
121 return;
122}
void ExchangeAllInPlace(vector< double > &xin, MPI_Comm const &icomm)
Definition vdop.cpp:107
void vdaxpy(std::vector< double > &x, std::vector< double > const &y, double alpha, std::vector< double > const &z)
Element-wise daxpy operation x(k) = y(k) + alpha*z(k).
Definition vdop.cpp:25
void vddiv(vector< double > &x, vector< double > const &y, vector< double > const &z)
Definition vdop.cpp:9
void ExchangeAll(vector< double > const &xin, vector< double > &yout, MPI_Comm const &icomm)
Definition vdop.cpp:87
double par_scalar(vector< double > const &x, vector< double > const &y, MPI_Comm const &icomm)
Definition vdop.cpp:77
double dscapr(std::vector< double > const &x, std::vector< double > const &y)
Calculates the Euclidean inner product of two vectors.
Definition vdop.cpp:40
bool CompareVectors(std::vector< double > const &x, int const n, double const y[], double const eps)
Compares an STL vector with POD vector.
Definition vdop.cpp:56