libROM  v1.0
Data-driven physical simulation library
Vector.cpp
1 
11 // Description: A simple, parallel Vector class with the utility needed to
12 // support the basis generation methods of this library. A
13 // distributed Vector has its rows distributed across processors.
14 
15 #include "Vector.h"
16 #include "utils/HDFDatabase.h"
17 
18 #include "mpi.h"
19 
20 #include <cmath>
21 #include <string.h>
22 
23 namespace CAROM {
24 
25 Vector::Vector() :
26  d_vec(NULL),
27  d_alloc_size(0),
28  d_distributed(false),
29  d_owns_data(true)
30 {}
31 
32 Vector::Vector(
33  int dim,
34  bool distributed) :
35  d_vec(NULL),
36  d_alloc_size(0),
37  d_distributed(distributed),
38  d_owns_data(true)
39 {
40  CAROM_VERIFY(dim > 0);
41  setSize(dim);
42  int mpi_init;
43  MPI_Initialized(&mpi_init);
44  if (mpi_init) {
45  MPI_Comm_size(MPI_COMM_WORLD, &d_num_procs);
46  }
47  else {
48  d_num_procs = 1;
49  }
50 }
51 
52 Vector::Vector(
53  double* vec,
54  int dim,
55  bool distributed,
56  bool copy_data) :
57  d_vec(NULL),
58  d_alloc_size(0),
59  d_distributed(distributed),
60  d_owns_data(copy_data)
61 {
62  CAROM_VERIFY(vec != 0);
63  CAROM_VERIFY(dim > 0);
64  if (copy_data) {
65  setSize(dim);
66  memcpy(d_vec, vec, d_alloc_size*sizeof(double));
67  }
68  else {
69  d_vec = vec;
70  d_alloc_size = dim;
71  d_dim = dim;
72  }
73  int mpi_init;
74  MPI_Initialized(&mpi_init);
75  if (mpi_init) {
76  MPI_Comm_size(MPI_COMM_WORLD, &d_num_procs);
77  }
78  else {
79  d_num_procs = 1;
80  }
81 }
82 
83 Vector::Vector(
84  const Vector& other) :
85  d_vec(NULL),
86  d_alloc_size(0),
87  d_distributed(other.d_distributed),
88  d_owns_data(true)
89 {
90  setSize(other.d_dim);
91  int mpi_init;
92  MPI_Initialized(&mpi_init);
93  if (mpi_init) {
94  MPI_Comm_size(MPI_COMM_WORLD, &d_num_procs);
95  }
96  else {
97  d_num_procs = 1;
98  }
99  memcpy(d_vec, other.d_vec, d_alloc_size*sizeof(double));
100 }
101 
103 {
104  if (d_owns_data && d_vec) {
105  delete [] d_vec;
106  }
107 }
108 
109 Vector&
111  const Vector& rhs)
112 {
113  d_distributed = rhs.d_distributed;
114  d_num_procs = rhs.d_num_procs;
115  setSize(rhs.d_dim);
116  memcpy(d_vec, rhs.d_vec, d_dim*sizeof(double));
117  return *this;
118 }
119 
120 Vector&
122  const Vector& rhs)
123 {
124  CAROM_VERIFY(d_dim == rhs.d_dim);
125  for(int i=0; i<d_dim; ++i) d_vec[i] += rhs.d_vec[i];
126  return *this;
127 }
128 
129 Vector&
131  const Vector& rhs)
132 {
133  CAROM_VERIFY(d_dim == rhs.d_dim);
134  for(int i=0; i<d_dim; ++i) d_vec[i] -= rhs.d_vec[i];
135  return *this;
136 }
137 
138 Vector&
139 Vector::operator = (const double& a)
140 {
141  for(int i=0; i<d_dim; ++i) d_vec[i] = a;
142  return *this;
143 }
144 
145 Vector&
146 Vector::operator *= (const double& a)
147 {
148  for(int i=0; i<d_dim; ++i) d_vec[i] *= a;
149  return *this;
150 }
151 
152 Vector&
153 Vector::transform(std::function<void(const int size, double* vector)>
154  transformer) {
155  transformer(d_dim, d_vec);
156  return *this;
157 }
158 
159 void
161  std::function<void(const int size, double* vector)> transformer) const {
162  result.setSize(d_dim);
163  result.d_distributed = d_distributed;
164  transformer(d_dim, result.d_vec);
165 }
166 
167 void
169  std::function<void(const int size, double* vector)> transformer) const {
170  // If the result has not been allocated then do so. Otherwise size it
171  // correctly.
172  if (result == 0) {
173  result = new Vector(d_dim, d_distributed);
174  }
175  else {
176  result->setSize(d_dim);
177  result->d_distributed = d_distributed;
178  }
179  transformer(d_dim, result->d_vec);
180 }
181 
182 Vector&
184  std::function<void(const int size, double* origVector, double* resultVector)>
185  transformer) {
186  Vector* origVector = new Vector(*this);
187  transformer(d_dim, origVector->d_vec, d_vec);
188  delete origVector;
189  return *this;
190 }
191 
192 void
194  std::function<void(const int size, double* origVector, double* resultVector)>
195  transformer) const {
196  Vector* origVector = new Vector(*this);
197  result.setSize(d_dim);
198  result.d_distributed = d_distributed;
199  transformer(d_dim, origVector->d_vec, result.d_vec);
200  delete origVector;
201 }
202 
203 void
205  std::function<void(const int size, double* origVector, double* resultVector)>
206  transformer) const {
207  // If the result has not been allocated then do so. Otherwise size it
208  // correctly.
209  Vector* origVector = new Vector(*this);
210  if (result == 0) {
211  result = new Vector(d_dim, d_distributed);
212  }
213  else {
214  result->setSize(d_dim);
215  result->d_distributed = d_distributed;
216  }
217  transformer(d_dim, origVector->d_vec, result->d_vec);
218  delete origVector;
219 }
220 
221 double
223  const Vector& other) const
224 {
225  CAROM_VERIFY(dim() == other.dim());
226  CAROM_ASSERT(distributed() == other.distributed());
227  double ip;
228  double local_ip = 0.0;
229  for (int i = 0; i < d_dim; ++i) {
230  local_ip += d_vec[i]*other.d_vec[i];
231  }
232  if (d_num_procs > 1 && d_distributed) {
233  MPI_Allreduce(&local_ip, &ip, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
234  }
235  else {
236  ip = local_ip;
237  }
238  return ip;
239 }
240 
241 double
243 {
244  double norm = sqrt(norm2());
245  return norm;
246 }
247 
248 double
250 {
251  double norm2 = inner_product(this);
252  return norm2;
253 }
254 
255 double
257 {
258  double Norm = norm();
259  for (int i = 0; i < d_dim; ++i) {
260  d_vec[i] /= Norm;
261  }
262  return Norm;
263 }
264 
265 void
267  const Vector& other,
268  Vector*& result) const
269 {
270  CAROM_ASSERT(result == 0 || result->distributed() == distributed());
271  CAROM_ASSERT(distributed() == other.distributed());
272  CAROM_VERIFY(dim() == other.dim());
273 
274  // If the result has not been allocated then do so. Otherwise size it
275  // correctly.
276  if (result == 0) {
277  result = new Vector(d_dim, d_distributed);
278  }
279  else {
280  result->setSize(d_dim);
281  }
282 
283  // Do the addition.
284  for (int i = 0; i < d_dim; ++i) {
285  result->d_vec[i] = d_vec[i] + other.d_vec[i];
286  }
287 }
288 
289 void
291  const Vector& other,
292  Vector& result) const
293 {
294  CAROM_ASSERT(result.distributed() == distributed());
295  CAROM_ASSERT(distributed() == other.distributed());
296  CAROM_VERIFY(dim() == other.dim());
297 
298  // Size result correctly.
299  result.setSize(d_dim);
300 
301  // Do the addition.
302  for (int i = 0; i < d_dim; ++i) {
303  result.d_vec[i] = d_vec[i] + other.d_vec[i];
304  }
305 }
306 
307 void
309  double factor,
310  const Vector& other,
311  Vector*& result) const
312 {
313  CAROM_ASSERT(result == 0 || result->distributed() == distributed());
314  CAROM_ASSERT(distributed() == other.distributed());
315  CAROM_VERIFY(dim() == other.dim());
316 
317  // If the result has not been allocated then do so. Otherwise size it
318  // correctly.
319  if (result == 0) {
320  result = new Vector(d_dim, d_distributed);
321  }
322  else {
323  result->setSize(d_dim);
324  }
325 
326  // Do the addition.
327  for (int i = 0; i < d_dim; ++i) {
328  result->d_vec[i] = d_vec[i] + factor*other.d_vec[i];
329  }
330 }
331 
332 void
334  double factor,
335  const Vector& other,
336  Vector& result) const
337 {
338  CAROM_ASSERT(result.distributed() == distributed());
339  CAROM_ASSERT(distributed() == other.distributed());
340  CAROM_VERIFY(dim() == other.dim());
341 
342  // Size result correctly.
343  result.setSize(d_dim);
344 
345  // Do the addition.
346  for (int i = 0; i < d_dim; ++i) {
347  result.d_vec[i] = d_vec[i] + factor*other.d_vec[i];
348  }
349 }
350 
351 void
353  double factor,
354  const Vector& other)
355 {
356  CAROM_ASSERT(distributed() == other.distributed());
357  CAROM_VERIFY(dim() == other.dim());
358 
359  // Do the addition.
360  for (int i = 0; i < d_dim; ++i) {
361  d_vec[i] += factor*other.d_vec[i];
362  }
363 }
364 
365 void
367  const Vector& other,
368  Vector*& result) const
369 {
370  CAROM_ASSERT(result == 0 || result->distributed() == distributed());
371  CAROM_ASSERT(distributed() == other.distributed());
372  CAROM_VERIFY(dim() == other.dim());
373 
374  // If the result has not been allocated then do so. Otherwise size it
375  // correctly.
376  if (result == 0) {
377  result = new Vector(d_dim, d_distributed);
378  }
379  else {
380  result->setSize(d_dim);
381  }
382 
383  // Do the subtraction.
384  for (int i = 0; i < d_dim; ++i) {
385  result->d_vec[i] = d_vec[i] - other.d_vec[i];
386  }
387 }
388 
389 void
391  const Vector& other,
392  Vector& result) const
393 {
394  CAROM_ASSERT(result.distributed() == distributed());
395  CAROM_ASSERT(distributed() == other.distributed());
396  CAROM_VERIFY(dim() == other.dim());
397 
398  // Size result correctly.
399  result.setSize(d_dim);
400 
401  // Do the subtraction.
402  for (int i = 0; i < d_dim; ++i) {
403  result.d_vec[i] = d_vec[i] - other.d_vec[i];
404  }
405 }
406 
407 void
409  double factor,
410  Vector*& result) const
411 {
412  CAROM_ASSERT(result == 0 || result->distributed() == distributed());
413 
414  // If the result has not been allocated then do so. Otherwise size it
415  // correctly.
416  if (result == 0) {
417  result = new Vector(d_dim, d_distributed);
418  }
419  else {
420  result->setSize(d_dim);
421  }
422 
423  // Do the multiplication.
424  for (int i = 0; i < d_dim; ++i) {
425  result->d_vec[i] = factor*d_vec[i];
426  }
427 }
428 
429 void
431  double factor,
432  Vector& result) const
433 {
434  CAROM_ASSERT(result.distributed() == distributed());
435 
436  // Size result correctly.
437  result.setSize(d_dim);
438 
439  // Do the multiplication.
440  for (int i = 0; i < d_dim; ++i) {
441  result.d_vec[i] = factor*d_vec[i];
442  }
443 }
444 
445 void
446 Vector::write(const std::string& base_file_name)
447 {
448  CAROM_ASSERT(!base_file_name.empty());
449 
450  int mpi_init;
451  MPI_Initialized(&mpi_init);
452  int rank;
453  if (mpi_init) {
454  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
455  }
456  else {
457  rank = 0;
458  }
459 
460  char tmp[100];
461  sprintf(tmp, ".%06d", rank);
462  std::string full_file_name = base_file_name + tmp;
463  HDFDatabase database;
464  database.create(full_file_name);
465 
466  sprintf(tmp, "distributed");
467  database.putInteger(tmp, d_distributed);
468  sprintf(tmp, "dim");
469  database.putInteger(tmp, d_dim);
470  sprintf(tmp, "data");
471  database.putDoubleArray(tmp, d_vec, d_dim);
472  database.close();
473 }
474 
475 void
476 Vector::print(const char * prefix) const
477 {
478  int my_rank;
479  const bool success = MPI_Comm_rank(MPI_COMM_WORLD, &my_rank);
480  CAROM_ASSERT(success);
481 
482  std::string filename_str = prefix + std::to_string(my_rank);
483  const char * filename = filename_str.c_str();
484  FILE * pFile = fopen(filename,"w");
485  for (int k = 0; k < d_dim; ++k) {
486  fprintf(pFile, " %25.20e\n", d_vec[k]);
487  }
488  fclose(pFile);
489 }
490 
491 void
492 Vector::read(const std::string& base_file_name)
493 {
494  CAROM_ASSERT(!base_file_name.empty());
495 
496  int mpi_init;
497  MPI_Initialized(&mpi_init);
498  int rank;
499  if (mpi_init) {
500  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
501  }
502  else {
503  rank = 0;
504  }
505 
506  char tmp[100];
507  sprintf(tmp, ".%06d", rank);
508  std::string full_file_name = base_file_name + tmp;
509  HDFDatabase database;
510  database.open(full_file_name, "r");
511 
512  sprintf(tmp, "distributed");
513  int distributed;
514  database.getInteger(tmp, distributed);
515  d_distributed = bool(distributed);
516  int dim;
517  sprintf(tmp, "dim");
518  database.getInteger(tmp, dim);
519  setSize(dim);
520  sprintf(tmp, "data");
521  database.getDoubleArray(tmp, d_vec, d_alloc_size);
522  d_owns_data = true;
523  if (mpi_init) {
524  MPI_Comm_size(MPI_COMM_WORLD, &d_num_procs);
525  }
526  else {
527  d_num_procs = 1;
528  }
529  database.close();
530 }
531 
532 void
533 Vector::local_read(const std::string& base_file_name, int rank)
534 {
535  CAROM_ASSERT(!base_file_name.empty());
536 
537  int mpi_init;
538  MPI_Initialized(&mpi_init);
539 
540  char tmp[100];
541  sprintf(tmp, ".%06d", rank);
542  std::string full_file_name = base_file_name + tmp;
543  HDFDatabase database;
544  database.open(full_file_name, "r");
545 
546  sprintf(tmp, "distributed");
547  int distributed;
548  database.getInteger(tmp, distributed);
549  d_distributed = bool(distributed);
550  int dim;
551  sprintf(tmp, "dim");
552  database.getInteger(tmp, dim);
553  setSize(dim);
554  sprintf(tmp, "data");
555  database.getDoubleArray(tmp, d_vec, d_alloc_size);
556  d_owns_data = true;
557  if (mpi_init) {
558  MPI_Comm_size(MPI_COMM_WORLD, &d_num_procs);
559  }
560  else {
561  d_num_procs = 1;
562  }
563  database.close();
564 }
565 
566 double Vector::localMin(int nmax)
567 {
568  const int n = nmax > 0 ? nmax : d_dim;
569  double v = d_vec[0];
570 
571  for (int i=1; i<n; ++i)
572  {
573  if (d_vec[i] < v)
574  v = d_vec[i];
575  }
576 
577  return v;
578 }
579 
580 int getCenterPoint(std::vector<Vector*>& points,
581  bool use_centroid)
582 {
583  int center_point;
584 
585  // get center-most point
586  // obtain the centroid and find the point closest to centroid
587  if (use_centroid)
588  {
589  Vector centroid(*points[0]);
590  for (int i = 1; i < points.size(); i++) {
591  centroid += *points[i];
592  }
593  centroid.mult(1.0 / points.size(), centroid);
594 
595  double min_dist_to_centroid;
596  for (int i = 0; i < points.size(); i++) {
597  Vector diff;
598  centroid.minus(*points[i], diff);
599  double dist_to_centroid = diff.norm();
600  if (i == 0 || dist_to_centroid < min_dist_to_centroid)
601  {
602  min_dist_to_centroid = dist_to_centroid;
603  center_point = i;
604  }
605  }
606  }
607 
608  // otherwise, do a brute force search to obtain
609  // the point closest to all other points
610  else
611  {
612  std::vector<double> dist_to_points(points.size(), 0);
613  for (int i = 0; i < points.size(); i++) {
614  for (int j = i + 1; j < points.size(); j++) {
615  Vector diff;
616  points[i]->minus(*points[j], diff);
617  double dist = diff.norm();
618  dist_to_points[i] += dist;
619  dist_to_points[j] += dist;
620  }
621  }
622 
623  double closest_dist_to_points = INT_MAX;
624  for (int i = 0; i < dist_to_points.size(); i++) {
625  if (dist_to_points[i] < closest_dist_to_points)
626  {
627  closest_dist_to_points = dist_to_points[i];
628  center_point = i;
629  }
630  }
631  }
632 
633  return center_point;
634 }
635 
636 int getCenterPoint(std::vector<Vector>& points,
637  bool use_centroid)
638 {
639  std::vector<Vector*> temp_points;
640  for (int i = 0; i < points.size(); i++)
641  {
642  temp_points.push_back(&points[i]);
643  }
644  return getCenterPoint(temp_points, use_centroid);
645 }
646 
647 int getClosestPoint(std::vector<Vector*>& points,
648  Vector* test_point)
649 {
650  int closest_point = 0;
651  double closest_dist_to_test_point = INT_MAX;
652  for (int i = 0; i < points.size(); i++) {
653  Vector diff;
654  test_point->minus(*points[i], diff);
655  double dist = diff.norm();
656  if (dist < closest_dist_to_test_point)
657  {
658  closest_dist_to_test_point = dist;
659  closest_point = i;
660  }
661  }
662 
663  return closest_point;
664 }
665 
666 int getClosestPoint(std::vector<Vector>& points,
667  Vector test_point)
668 {
669  std::vector<Vector*> temp_points;
670  for (int i = 0; i < points.size(); i++)
671  {
672  temp_points.push_back(&points[i]);
673  }
674  return getClosestPoint(temp_points, &test_point);
675 }
676 
677 }
void getInteger(const std::string &key, int &data)
Reads an integer associated with the supplied key from the currently open database file.
Definition: Database.h:182
void putInteger(const std::string &key, int data)
Writes an integer associated with the supplied key to currently open database file.
Definition: Database.h:95
virtual bool create(const std::string &file_name, const MPI_Comm comm=MPI_COMM_NULL) override
Creates a new HDF5 database file with the supplied name.
Definition: HDFDatabase.cpp:45
virtual bool open(const std::string &file_name, const std::string &type, const MPI_Comm comm=MPI_COMM_NULL) override
Opens an existing HDF5 database file with the supplied name.
Definition: HDFDatabase.cpp:76
virtual void getDoubleArray(const std::string &key, double *data, int nelements, const bool distributed=false)
Reads an array of doubles associated with the supplied key from the currently open HDF5 database file...
virtual void putDoubleArray(const std::string &key, const double *const data, int nelements, const bool distributed=false)
Writes an array of doubles associated with the supplied key to the currently open HDF5 database file.
virtual bool close()
Closes the currently open HDF5 database file.
void setSize(int dim)
Sets the length of the vector and reallocates storage if needed. All values are initialized to zero.
Definition: Vector.h:231
Vector * mult(double factor) const
Multiplies this by the supplied constant and returns the result.
Definition: Vector.h:610
void read(const std::string &base_file_name)
read Vector from (a) HDF file(s).
Definition: Vector.cpp:492
Vector & operator-=(const Vector &rhs)
Subtraction operator.
Definition: Vector.cpp:130
bool distributed() const
Returns true if the Vector is distributed.
Definition: Vector.h:255
void print(const char *prefix) const
print Vector into (a) ascii file(s).
Definition: Vector.cpp:476
double norm() const
Form the norm of this.
Definition: Vector.cpp:242
Vector * minus(const Vector &other) const
Subtracts other and this and returns the result, reference version.
Definition: Vector.h:538
double normalize()
Normalizes the Vector and returns its norm.
Definition: Vector.cpp:256
void write(const std::string &base_file_name)
write Vector into (a) HDF file(s).
Definition: Vector.cpp:446
Vector * plus(const Vector &other) const
Adds other and this and returns the result, reference version.
Definition: Vector.h:349
Vector & operator=(const Vector &rhs)
Assignment operator.
Definition: Vector.cpp:110
double inner_product(const Vector &other) const
Inner product, reference form.
Definition: Vector.cpp:222
~Vector()
Destructor.
Definition: Vector.cpp:102
void plusEqAx(double factor, const Vector &other)
Adds factor*other to this, reference version.
Definition: Vector.cpp:352
double norm2() const
Form the squared norm of this.
Definition: Vector.cpp:249
Vector * plusAx(double factor, const Vector &other)
Adds factor*other and this and returns the result, reference version.
Definition: Vector.h:424
Vector & transform(std::function< void(const int size, double *vector)> transformer)
Transform the vector using a supplied function.
Definition: Vector.cpp:153
Vector & operator+=(const Vector &rhs)
Addition operator.
Definition: Vector.cpp:121
double localMin(int nmax=0)
Compute the local minimum of this.
Definition: Vector.cpp:566
void local_read(const std::string &base_file_name, int rank)
read read a single rank of a distributed Vector from (a) HDF file(s).
Definition: Vector.cpp:533
int dim() const
Returns the dimension of the Vector on this processor.
Definition: Vector.h:266
Vector & operator*=(const double &a)
Scaling operator.
Definition: Vector.cpp:146
int getClosestPoint(std::vector< Vector * > &points, Vector *test_point)
Get closest point to a test point among a group of points.
Definition: Vector.cpp:647
int getCenterPoint(std::vector< Vector * > &points, bool use_centroid)
Get center point of a group of points.
Definition: Vector.cpp:580