OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
floatarray.C
Go to the documentation of this file.
1 /*
2  *
3  * ##### ##### ###### ###### ### ###
4  * ## ## ## ## ## ## ## ### ##
5  * ## ## ## ## #### #### ## # ##
6  * ## ## ## ## ## ## ## ##
7  * ## ## ## ## ## ## ## ##
8  * ##### ##### ## ###### ## ##
9  *
10  *
11  * OOFEM : Object Oriented Finite Element Code
12  *
13  * Copyright (C) 1993 - 2013 Borek Patzak
14  *
15  *
16  *
17  * Czech Technical University, Faculty of Civil Engineering,
18  * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19  *
20  * This library is free software; you can redistribute it and/or
21  * modify it under the terms of the GNU Lesser General Public
22  * License as published by the Free Software Foundation; either
23  * version 2.1 of the License, or (at your option) any later version.
24  *
25  * This program is distributed in the hope that it will be useful,
26  * but WITHOUT ANY WARRANTY; without even the implied warranty of
27  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28  * Lesser General Public License for more details.
29  *
30  * You should have received a copy of the GNU Lesser General Public
31  * License along with this library; if not, write to the Free Software
32  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33  */
34 
35 #include "floatarray.h"
36 #include "intarray.h"
37 #include "floatmatrix.h"
38 #include "mathfem.h"
39 #include "error.h"
40 #include "datastream.h"
41 #include "mathfem.h"
42 
43 #include <cstdarg>
44 #include <cstdlib>
45 #include <cstring>
46 #include <ostream>
47 #include <memory>
48 #include <numeric>
49 #include <cmath>
50 #include <iostream>
51 #include <fstream>
52 #include <iomanip>
53 
54 #define FAST_RESIZE(newsize) \
55  if ( (newsize) < this->giveSize() ) { \
56  this->values.resize((newsize)); \
57  } else if ( (newsize) > this->giveSize() ) { \
58  this->values.assign((newsize), 0.); \
59  }
60 
61 #ifdef __LAPACK_MODULE
62 extern "C" {
63  extern void dgemv_(const char *trans, const int *m, const int *n, const double *alpha, const double *a, const int *lda, const double *x,
64  const int *incx, const double *beta, double *y, const int *incy, int aColumns, int xSize, int ySize);
65  // Y = Y + alpha * X
66  extern void daxpy_(const int *n, const double *alpha, const double *x, const int *incx, double *y, const int *incy, int xsize, int ysize);
67 }
68 #endif
69 
70 namespace oofem {
71 
73 {
74  for(double val : values) {
75  if(!std::isfinite(val)) {
76  return false;
77  }
78  }
79 
80  return true;
81 }
82 
83 #ifdef DEBUG
84 double &
86 {
87  if ( i >= this->giveSize() ) {
88  OOFEM_ERROR("array error on index : %d >= %d", i, this->giveSize());
89  }
90  return values [ i ];
91 }
92 
93 const double &
94 FloatArray :: operator() (int i) const
95 {
96  if ( i >= this->giveSize() ) {
97  OOFEM_ERROR("array error on index : %d >= %d", i, this->giveSize());
98  }
99  return values [ i ];
100 }
101 double &
103 {
104  if ( i >= this->giveSize() ) {
105  OOFEM_ERROR("array error on index : %d >= %d", i, this->giveSize());
106  }
107  return values [ i ];
108 }
109 
110 const double &
111 FloatArray :: operator[] (int i) const
112 {
113  if ( i >= this->giveSize() ) {
114  OOFEM_ERROR("array error on index : %d >= %d", i, this->giveSize());
115  }
116  return values [ i ];
117 }
118 
119 double &FloatArray :: at(int i)
120 {
121  this->checkBounds(i);
122  return values [ i - 1 ];
123 }
124 
125 double FloatArray :: at(int i) const
126 {
127  this->checkBounds(i);
128  return values [ i - 1 ];
129 }
130 
131 #endif
132 void FloatArray :: checkBounds(int i) const
133 // Checks that the receiver's size is not smaller than 'i'.
134 {
135  if ( i <= 0 ) {
136  OOFEM_ERROR("array error on index : %d <= 0", i);
137  }
138 
139  if ( i > this->giveSize() ) {
140  OOFEM_ERROR("array error on index : %d > %d", i, this->giveSize());
141  }
142 }
143 
144 
145 void
147 {
148  FAST_RESIZE(b.giveSize());
149 
150  for ( int i = 0; i < this->giveSize(); ++i ) {
151  (*this) [ i ] = s * b [ i ];
152  }
153 }
154 
155 
157 // Performs the operation a=a+b, where a stands for the receiver. If the
158 // receiver's size is 0, adjusts its size to that of b. Returns the
159 // receiver.
160 {
161  if ( b.isEmpty() ) {
162  return;
163  }
164 
165  if ( !this->giveSize() ) {
166  * this = b;
167  return;
168  }
169 
170 # ifdef DEBUG
171  if ( this->giveSize() != b.giveSize() ) {
172  OOFEM_ERROR("dimension mismatch in a[%d]->add(b[%d])", this->giveSize(), b.giveSize());
173  }
174 
175 # endif
176 
177 #ifdef __LAPACK_MODULE
178  int inc = 1;
179  double s = 1.;
180  int size = this->giveSize();
181  daxpy_(& size, & s, b.givePointer(), & inc, this->givePointer(), & inc, b.giveSize(), this->giveSize());
182 #else
183  for ( int i = 0; i < this->giveSize(); i++ ) {
184  (*this) [ i ] += b [ i ];
185  }
186 #endif
187 }
188 
189 
190 void FloatArray :: add(double offset)
191 {
192  for ( double &x: *this ) {
193  x += offset;
194  }
195 }
196 
197 
198 void FloatArray :: add(double factor, const FloatArray &b)
199 // Performs the operation a=a+factor*b, where a stands for the receiver. If the
200 // receiver's size is 0, adjusts its size to that of b.
201 {
202  if ( this->isEmpty() ) {
203  this->beScaled(factor, b);
204  return;
205  }
206 
207 # ifdef DEBUG
208  if ( this->giveSize() != b.giveSize() ) {
209  OOFEM_ERROR("dimension mismatch in a[%d]->add(b[%d])", this->giveSize(), b.giveSize());
210  }
211 
212 # endif
213 
214 #ifdef __LAPACK_MODULE
215  int inc = 1;
216  int size = this->giveSize();
217  daxpy_(& size, & factor, b.givePointer(), & inc, this->givePointer(), & inc, b.giveSize(), this->giveSize());
218 #else
219  for ( int i = 0; i < this->giveSize(); ++i ) {
220  (*this) [ i ] += factor * b [ i ];
221  }
222 #endif
223 }
224 
225 
226 void FloatArray :: plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
227 // Performs the operation a += b^T . s * dV
228 {
229  int nRows = b.giveNumberOfRows();
230  int nColumns = b.giveNumberOfColumns();
231 
232  if ( this->isEmpty() ) {
233  this->values.assign( nColumns, 0. );
234  }
235 
236 # ifdef DEBUG
237  if ( this->giveSize() != b.giveNumberOfColumns() ) {
238  OOFEM_ERROR( "dimension mismatch in a[%d] and b[%d, *]", this->giveSize(), b.giveNumberOfColumns() );
239  }
240 # endif
241 
242 #ifdef __LAPACK_MODULE
243  double beta = 1.;
244  int inc = 1;
245  dgemv_("t", & nRows, & nColumns, & dV, b.givePointer(), & nRows, s.givePointer(), & inc, & beta, this->givePointer(), & inc, nColumns, nColumns, nRows);
246 #else
247  for ( int i = 1; i <= nColumns; i++ ) {
248  double sum = 0.;
249  for ( int j = 1; j <= nRows; j++ ) {
250  sum += b.at(j, i) * s.at(j);
251  }
252  this->at(i) += sum * dV;
253  }
254 #endif
255 }
256 
257 
259 // Performs the operation a=a-src, where a stands for the receiver. If the
260 // receiver's size is 0, adjusts its size to that of src.
261 {
262  if ( src.isEmpty() ) {
263  return;
264  }
265 
266  if ( this->isEmpty() ) {
267  FAST_RESIZE(src.giveSize());
268  for ( int i = 0; i < this->giveSize(); ++i ) {
269  (*this) [ i ] = -src [ i ];
270  }
271 
272  return;
273  }
274 
275 # ifdef DEBUG
276  if ( this->giveSize() != src.giveSize() ) {
277  OOFEM_ERROR("dimension mismatch in a[%d]->add(b[%d])", this->giveSize(), src.giveSize());
278  }
279 
280 # endif
281 
282  for ( int i = 0; i < this->giveSize(); ++i ) {
283  (*this) [ i ] -= src [ i ];
284  }
285 }
286 
287 
289 {
290  int n = a.giveSize();
291 
292  if ( a.giveSize() == 0 ) {
293  *this = b;
294  return;
295  } else if ( b.giveSize() == 0 ) {
296  *this = a;
297  return;
298  }
299 
300 # ifdef DEBUG
301  if ( n != b.giveSize() ) {
302  OOFEM_ERROR("dimension mismatch in beMaxOf(a[%d],b[%d])", n, b.giveSize());
303  }
304 
305 # endif
306 
307  FAST_RESIZE(n);
308 
309  for ( int i = 0; i < n; i++ ) {
310  (*this) [ i ] = max( a [ i ], b [ i ] );
311  }
312 }
313 
314 
316 {
317  int n = a.giveSize();
318 
319  if ( a.giveSize() == 0 ) {
320  *this = b;
321  return;
322  } else if ( b.giveSize() == 0 ) {
323  *this = a;
324  return;
325  }
326 
327 # ifdef DEBUG
328  if ( n != b.giveSize() ) {
329  OOFEM_ERROR("dimension mismatch in beMinOf(a[%d],b[%d])", n, b.giveSize());
330  }
331 
332 # endif
333 
334  FAST_RESIZE(n);
335  for ( int i = 0; i < n; i++ ) {
336  (*this) [ i ] = min( a [ i ], b [ i ] );
337  }
338 }
339 
340 
342 {
343 #ifdef DEBUG
344  if ( a.giveSize() != b.giveSize() ) {
345  OOFEM_ERROR("size mismatch (%d : %d)", a.giveSize(), b.giveSize());
346  }
347 
348 #endif
349 #if 0
350  FAST_RESIZE(a.giveSize());
351  for ( int i = 0; i < this->giveSize(); ++i ) {
352  (*this) [ i ] = a [ i ] - b [ i ];
353  }
354 #else
355  this->values.reserve(a.giveSize());
356  this->values.resize(0);
357  for ( int i = 0; i < a.giveSize(); ++i ) {
358  this->values.push_back( a[i] - b[i] );
359  }
360 
361 #endif
362 }
363 
364 void FloatArray :: beDifferenceOf(const FloatArray &a, const FloatArray &b, int n)
365 {
366 #ifdef DEBUG
367  if ( a.giveSize() < n || b.giveSize() < n ) {
368  OOFEM_ERROR("wrong size ", a.giveSize(), b.giveSize());
369  }
370 
371 #endif
372  FAST_RESIZE(n);
373  for ( int i = 0; i < n; ++i ) {
374  (*this) [ i ] = a [ i ] - b [ i ];
375  }
376 }
377 
378 
379 void FloatArray :: beSubArrayOf(const FloatArray &src, const IntArray &indx)
380 //
381 // Returns a subVector of the receiver constructed from an index array.
382 // this->at(i) = src.at(indx->at(i))
383 //
384 {
385 #ifdef DEBUG
386  if ( indx.maximum() > src.giveSize() || indx.minimum() < 1 ) {
387  OOFEM_ERROR("index points outside of source");
388  }
389 #endif
390 
391  int n = indx.giveSize();
392  FAST_RESIZE(n);
393  for ( int i = 1; i <= n; i++ ) {
394  this->at(i) = src.at( indx.at(i) );
395  }
396 }
397 
398 
399 void FloatArray :: addSubVector(const FloatArray &src, int si)
400 {
401  int reqSize, n = src.giveSize();
402 
403  si--;
404  reqSize = si + n;
405  if ( this->giveSize() < reqSize ) {
406  this->resizeWithValues(reqSize);
407  }
408 
409  for ( int i = 0; i < n; i++ ) {
410  (*this) [si + i] += src [ i ];
411  }
412 }
413 
414 
416 {
417 # if DEBUG
418  // check proper bounds
419  if ( ( v1.giveSize() != 3 ) || ( v2.giveSize() != 3 ) ) {
420  OOFEM_ERROR("size mismatch, size is not equal to 3");
421  }
422 # endif
423 
424  FAST_RESIZE(3);
425 
426  this->at(1) = v1.at(2) * v2.at(3) - v1.at(3) * v2.at(2);
427  this->at(2) = v1.at(3) * v2.at(1) - v1.at(1) * v2.at(3);
428  this->at(3) = v1.at(1) * v2.at(2) - v1.at(2) * v2.at(1);
429 }
430 
432 {
433  int index = 1;
434  if ( !this->giveSize() ) {
435  return -1;
436  }
437  double val = (*this) [ 0 ];
438  for ( int i = 1; i < this->giveSize(); i++ ) {
439  if ( val > (*this) [ i ] ) {
440  val = (*this) [ i ];
441  index = i + 1;
442  }
443  }
444  return index;
445 }
446 
448 {
449  int index = 1;
450  if ( !this->giveSize() ) {
451  return -1;
452  }
453  double val = (*this) [ 0 ];
454  for ( int i = 1; i < this->giveSize(); i++ ) {
455  if ( val < (*this) [ i ] ) {
456  val = (*this) [ i ];
457  index = i + 1;
458  }
459  }
460  return index;
461 }
462 
463 double FloatArray :: dotProduct(const FloatArray &x) const
464 {
465 # ifdef DEBUG
466  if ( this->giveSize() != x.giveSize() ) {
467  OOFEM_ERROR("dimension mismatch in a[%d]->dotProduct(b[%d])", this->giveSize(), x.giveSize());
468  }
469 
470 # endif
471 
472  return std::inner_product(this->begin(), this->end(), x.begin(), 0.);
473 }
474 
475 
476 double FloatArray :: dotProduct(const FloatArray &x, int size) const
477 {
478 # ifdef DEBUG
479  if ( size > this->giveSize() || size > x.giveSize() ) {
480  OOFEM_ERROR("dimension mismatch in a[%d]->dotProduct(b[%d])", this->giveSize(), x.giveSize());
481  }
482 
483 # endif
484 
485  return std::inner_product(this->begin(), this->begin()+size, x.begin(), 0.);
486 }
487 
488 
489 double FloatArray :: distance(const FloatArray &x) const
490 {
491  return sqrt( this->distance_square(x) );
492 }
493 
494 double FloatArray :: distance(const FloatArray &iP1, const FloatArray &iP2, double &oXi, double &oXiUnbounded) const
495 {
496  return sqrt( distance_square(iP1, iP2, oXi, oXiUnbounded) );
497 }
498 
499 double FloatArray :: distance_square(const FloatArray &iP1, const FloatArray &iP2, double &oXi, double &oXiUnbounded) const
500 {
501  const double l2 = iP1.distance_square(iP2);
502 
503  if ( l2 > 0.0 ) {
504 
505  const double s = (dot(*this, iP2) - dot(*this, iP1) ) - ( dot(iP1, iP2) - dot(iP1, iP1) );
506 
507  if ( s < 0.0 ) {
508  // X is closest to P1
509  oXi = 0.0;
510  oXiUnbounded = s/l2;
511  return this->distance_square(iP1);
512  } else {
513  if ( s > l2 ) {
514  // X is closest to P2
515  oXi = 1.0;
516  oXiUnbounded = s/l2;
517  return this->distance_square(iP2);
518  } else {
519  oXi = s / l2;
520  oXiUnbounded = s/l2;
521  const FloatArray q = ( 1.0 - oXi ) * iP1 + oXi * iP2;
522  return this->distance_square(q);
523  }
524  }
525  } else {
526  // If the points P1 and P2 coincide,
527  // we can compute the distance to any
528  // of these points.
529  oXi = 0.5;
530  oXiUnbounded = 0.5;
531  return this->distance_square(iP1);
532  }
533 }
534 
535 
536 double FloatArray :: distance_square(const FloatArray &from) const
537 // returns distance between receiver and from from
538 // computed using generalized pythagorean formulae
539 {
540  double dist = 0.;
541  int s = min(this->giveSize(), from.giveSize());
542  for ( int i = 1; i <= s; ++i ) {
543  double dx = this->at(i) - from.at(i);
544  dist += dx * dx;
545  }
546 
547  return dist;
548 }
549 
550 
551 void FloatArray :: assemble(const FloatArray &fe, const IntArray &loc)
552 // Assembles the array fe (typically, the load vector of a finite
553 // element) to the receiver, using loc as location array.
554 {
555  int n = fe.giveSize();
556 # ifdef DEBUG
557  if ( n != loc.giveSize() ) {
558  OOFEM_ERROR("dimensions of 'fe' (%d) and 'loc' (%d) mismatch", fe.giveSize(), loc.giveSize() );
559  }
560 
561 # endif
562 
563  for ( int i = 1; i <= n; i++ ) {
564  int ii = loc.at(i);
565  if ( ii ) { // if non 0 coefficient,
566  this->at(ii) += fe.at(i);
567  }
568  }
569 }
570 
571 
573 // Assembles the array fe (typically, the load vector of a finite
574 // element) to the receiver, using loc as location array.
575 {
576  int n = fe.giveSize();
577 # ifdef DEBUG
578  if ( n != loc.giveSize() ) {
579  OOFEM_ERROR("dimensions of 'fe' (%d) and 'loc' (%d) mismatch", fe.giveSize(), loc.giveSize() );
580  }
581 
582 # endif
583 
584  for ( int i = 1; i <= n; i++ ) {
585  int ii = loc.at(i);
586  if ( ii ) { // if non 0 coefficient,
587  this->at(ii) += fe.at(i) * fe.at(i);
588  }
589  }
590 }
591 
592 
594 // Expands the receiver if loc points to coefficients beyond the size of
595 // the receiver.
596 {
597  int high = 0;
598  for ( int p : loc ) {
599  high = max( high, p );
600  }
601 
602  if ( high > this->giveSize() ) { // receiver must be expanded
603  this->values.resize(high);
604  }
605 }
606 
607 
609 {
610  this->values.reserve(s);
611  this->values.clear();
612 }
613 
614 
615 void FloatArray :: resizeWithValues(int n, int allocChunk)
616 {
617 #ifdef DEBUG
618  if ( allocChunk < 0 ) {
619  OOFEM_FATAL("allocChunk must be non-negative; %d", allocChunk);
620  }
621 
622 #endif
623 
624  if ( allocChunk > 0 && (int)this->values.capacity() < n ) {
625  this->values.reserve(n + allocChunk);
626  }
627 
628  this->values.resize(n);
629 }
630 
632 {
633  this->values.resize(n);
635  //this->values.assign(n, 0.);
636 }
637 
638 
640 {
641  this->values.assign(n, 0.);
642  this->values.shrink_to_fit();
643 }
644 
645 
647 {
648  for ( double x : *this ) {
649  if ( x != 0. ) {
650  return false;
651  }
652  }
653 
654  return true;
655 }
656 
657 
659 {
660  std::fill(this->begin(), this->end(), 0.);
661 }
662 
663 
665 {
666  this->values.insert(this->end(), a.begin(), a.end());
667 }
668 
669 
670 void FloatArray :: append(double a)
671 {
672  this->values.push_back(a);
673 }
674 
675 
676 void FloatArray :: beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
677 // Stores the product of aMatrix * anArray in to receiver
678 {
679  int nColumns = aMatrix.giveNumberOfColumns();
680  int nRows = aMatrix.giveNumberOfRows();
681 
682  FAST_RESIZE(nRows);
683 
684 # ifdef DEBUG
685  if ( aMatrix.giveNumberOfColumns() != anArray.giveSize() ) {
686  OOFEM_ERROR("dimension mismatch (%d, %d) . (%d)",
687  aMatrix.giveNumberOfRows(), aMatrix.giveNumberOfColumns(), anArray.giveSize());
688  }
689 # endif
690 
691 #ifdef __LAPACK_MODULE
692  double alpha = 1., beta = 0.;
693  int inc = 1;
694  dgemv_("n", & nRows, & nColumns, & alpha, aMatrix.givePointer(), & nRows, anArray.givePointer(), & inc, & beta, this->givePointer(), & inc, nColumns, nColumns, nRows);
695 #else
696  for ( int i = 1; i <= nRows; i++ ) {
697  double sum = 0.;
698  for ( int j = 1; j <= nColumns; j++ ) {
699  sum += aMatrix.at(i, j) * anArray.at(j);
700  }
701 
702  this->at(i) = sum;
703  }
704 #endif
705 }
706 
707 
708 void FloatArray :: beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
709 // Stores the product of aMatrix^T * anArray in to receiver
710 {
711  int nRows = aMatrix.giveNumberOfRows();
712  int nColumns = aMatrix.giveNumberOfColumns();
713 
714 # ifdef DEBUG
715  if ( aMatrix.giveNumberOfRows() != anArray.giveSize() ) {
716  OOFEM_ERROR( "dimension mismatch, matrix rows = %d, array size = %d", aMatrix.giveNumberOfRows(), anArray.giveSize() );
717  }
718 
719 # endif
720  FAST_RESIZE(nColumns);
721 
722 #ifdef __LAPACK_MODULE
723  double alpha = 1., beta = 0.;
724  int inc = 1;
725  dgemv_("t", & nRows, & nColumns, & alpha, aMatrix.givePointer(), & nRows, anArray.givePointer(), & inc, & beta, this->givePointer(), & inc, nColumns, nColumns, nRows);
726 #else
727  for ( int i = 1; i <= nColumns; i++ ) {
728  double sum = 0.;
729  for ( int j = 1; j <= nRows; j++ ) {
730  sum += aMatrix.at(j, i) * anArray.at(j);
731  }
732 
733  this->at(i) = sum;
734  }
735 #endif
736 }
737 
738 
740 {
741  for ( double &x: *this ) {
742  x = -x;
743  }
744 }
745 
746 
748 // Prints the receiver on screen.
749 {
750  printf("FloatArray of size : %d \n", this->giveSize());
751  for ( double x: *this ) {
752  printf( "%10.3e ", x );
753  }
754 
755  printf("\n");
756 }
757 
758 
759 void FloatArray :: printYourself(const std::string &name) const
760 // Prints the receiver on screen.
761 {
762  printf("%s (%d): \n", name.c_str(), this->giveSize());
763  for ( double x: *this ) {
764  printf( "%10.3e ", x );
765  }
766 
767  printf("\n");
768 }
769 
770 void FloatArray :: printYourselfToFile(const std::string filename, const bool showDimensions) const
771 // Prints the receiver to file.
772 {
773  std :: ofstream arrayfile (filename);
774  if (arrayfile.is_open()) {
775  if (showDimensions)
776  arrayfile << "FloatArray of size : " << this->giveSize() << "\n";
777  arrayfile << std::scientific << std::right << std::setprecision(3);
778  for ( double x: *this ) {
779  arrayfile << std::setw(10) << x << "\t";
780  }
781  arrayfile.close();
782  } else {
783  OOFEM_ERROR("Failed to write to file");
784  }
785 }
786 
787 void FloatArray :: pY() const
788 // Prints the receiver on screen with higher accuracy than printYourself.
789 {
790  printf("[");
791  for ( double x: *this ) {
792  printf( "%20.14e; ", x );
793  }
794 
795  printf("];\n");
796 }
797 
798 
800 // Returns the receiver 'a' rotated according the change-of-base matrix r.
801 // If mode = 't', the method performs the operation a = r(transp) * a .
802 // If mode = 'n', the method performs the operation a = r * a .
803 {
804  FloatArray rta;
805 
806  if ( mode == 't' ) {
807  rta.beTProductOf(r, * this);
808  } else if ( mode == 'n' ) {
809  rta.beProductOf(r, * this);
810  } else {
811  OOFEM_ERROR("unsupported mode");
812  }
813 
814  * this = rta;
815 }
816 
817 
818 void FloatArray :: times(double factor)
819 // Multiplies every coefficient of the receiver by factor.
820 {
821  //std::transform(this->begin(), this->end(), this->begin(), std::bind2nd(std::multiplies<double>(), factor));
822  for ( double &x : *this ) {
823  x *= factor;
824  }
825 }
826 
827 
829 {
830  double norm = this->computeNorm();
831  if ( norm < 1.e-80 ) {
832  OOFEM_ERROR("cannot norm receiver, norm is too small");
833  }
834 
835  this->times(1. / norm);
836  return norm;
837 }
838 
839 
841 {
842  return sqrt( this->computeSquaredNorm() );
843 }
844 
845 
847 {
848  return std::inner_product(this->begin(), this->end(), this->begin(), 0.);
849 }
850 
851 
852 double FloatArray :: sum() const
853 {
854  return std::accumulate(this->begin(), this->end(), 0.);
855 }
856 
857 
858 double FloatArray :: product() const
859 {
860  return std::accumulate(this->begin(), this->end(), 1.0, [](double a, double b) { return a*b; });
861 }
862 
863 
864 void FloatArray :: copySubVector(const FloatArray &src, int si)
865 {
866  si--;
867  this->resizeWithValues(si + src.giveSize());
868  std :: copy( src.begin(), src.end(), this->begin() + si);
869 }
870 
871 
873 // writes receiver's binary image into stream
874 // use id to distinguish some instances
875 // return value >0 success
876 // =0 file i/o error
877 {
878  // write size
879  int size = this->giveSize();
880  if ( !stream.write(size) ) {
881  return CIO_IOERR;
882  }
883 
884  // write raw data
885  if ( size ) {
886  if ( !stream.write(this->givePointer(), size) ) {
887  return CIO_IOERR;
888  }
889  }
890 
891  // return result back
892  return CIO_OK;
893 }
894 
896 // reads receiver from stream
897 // warning - overwrites existing data!
898 // returns 0 if file i/o error
899 // -1 if id od class id is not correct
900 {
901  // read size
902  int size;
903  if ( !stream.read(size) ) {
904  return CIO_IOERR;
905  }
906 
907  this->values.resize(size);
908 
909  // read raw data
910  if ( size ) {
911  if ( !stream.read(this->givePointer(), size) ) {
912  return CIO_IOERR;
913  }
914  }
915 
916  // return result back
917  return CIO_OK;
918 }
919 
920 
922 {
923  return buff.givePackSizeOfInt(1) + buff.givePackSizeOfDouble(this->giveSize());
924 }
925 
926 // IML compat
927 
928 FloatArray &FloatArray :: operator = ( const double & val )
929 {
930  std::fill(this->begin(), this->begin(), val);
931  return * this;
932 }
933 
934 FloatArray &operator *= ( FloatArray & x, const double & a )
935 {
936  x.times(a);
937  return x;
938 }
939 
940 FloatArray operator *( const double & a, const FloatArray & x )
941 {
942  FloatArray result;
943  result.beScaled(a, x);
944  return result;
945 }
946 
947 FloatArray operator *( const FloatArray & x, const double & a )
948 {
949  FloatArray result;
950  result.beScaled(a, x);
951  return result;
952 }
953 
955 {
956  FloatArray result(x);
957  result.add(y);
958  return result;
959 }
960 
962 {
963  FloatArray result;
964  result.beDifferenceOf(x, y);
965  return result;
966 }
967 
969 {
970  x.add(y);
971  return x;
972 }
973 
975 {
976  x.subtract(y);
977  return x;
978 }
979 
980 double dot(const FloatArray &x, const FloatArray &y)
981 {
982  return x.dotProduct(y);
983 }
984 
985 double norm(const FloatArray &x)
986 {
987  return x.computeNorm();
988 }
989 
990 // End of IML compat
991 
993 {
994  // Rewrites the matrix on vector form, order: 11, 22, 33, 23, 13, 12, 32, 31, 21
995 # ifdef DEBUG
996  if ( aMatrix.giveNumberOfColumns() != 3 || aMatrix.giveNumberOfColumns() != 3 ) {
997  OOFEM_ERROR("matrix dimension is not 3x3");
998  }
999 
1000 # endif
1001  *this = {
1002  aMatrix.at(1, 1),
1003  aMatrix.at(2, 2),
1004  aMatrix.at(3, 3),
1005  aMatrix.at(2, 3),
1006  aMatrix.at(1, 3),
1007  aMatrix.at(1, 2),
1008  aMatrix.at(3, 2),
1009  aMatrix.at(3, 1),
1010  aMatrix.at(2, 1)
1011  };
1012 }
1013 
1015 {
1016  // Revrites a symmetric strain matrix on reduced vector form, order: 11, 22, 33, 23, 13, 12
1017  // shear components are multiplied with a factor 2
1018 # ifdef DEBUG
1019  if ( aMatrix.giveNumberOfColumns() != 3 || aMatrix.giveNumberOfColumns() != 3 ) {
1020  OOFEM_ERROR("matrix dimension is not 3x3");
1021  }
1022 # endif
1023 
1024  *this = {
1025  aMatrix.at(1, 1),
1026  aMatrix.at(2, 2),
1027  aMatrix.at(3, 3),
1028  ( aMatrix.at(2, 3) + aMatrix.at(3, 2) ),
1029  ( aMatrix.at(1, 3) + aMatrix.at(3, 1) ),
1030  ( aMatrix.at(1, 2) + aMatrix.at(2, 1) )
1031  };
1032 }
1033 
1034 
1036 {
1037  // Revrites the matrix on vector form (symmetrized matrix used), order: 11, 22, 33, 23, 13, 12
1038 # ifdef DEBUG
1039  if ( aMatrix.giveNumberOfColumns() != 3 || aMatrix.giveNumberOfColumns() != 3 ) {
1040  OOFEM_ERROR("matrix dimension is not 3x3");
1041  }
1042 
1043 # endif
1044 
1045  *this = {
1046  aMatrix.at(1, 1),
1047  aMatrix.at(2, 2),
1048  aMatrix.at(3, 3),
1049  0.5 * ( aMatrix.at(2, 3) + aMatrix.at(3, 2) ),
1050  0.5 * ( aMatrix.at(1, 3) + aMatrix.at(3, 1) ),
1051  0.5 * ( aMatrix.at(1, 2) + aMatrix.at(2, 1) )
1052  };
1053 }
1054 
1056 {
1057  // OOFEM: 11, 22, 33, 23, 13, 12, 32, 31, 21
1058  // UMAT: 11, 22, 33, 12, 13, 23, 32, 21, 31
1059 
1060  if ( this->giveSize() == 6 ) {
1061  std :: swap( this->at(4), this->at(6) );
1062  } else if ( this->giveSize() == 9 ) {
1063  // OOFEM: 11, 22, 33, 23, 13, 12, 32, 31, 21
1064  // UMAT: 11, 22, 33, 12, 13, 23, 32, 21, 31
1065  const int abq2oo [ 9 ] = {
1066  1, 2, 3, 6, 5, 4, 7, 9, 8
1067  };
1068 
1069  FloatArray tmp(9);
1070  for ( int i = 0; i < 9; i++ ) {
1071  tmp(i) = this->at(abq2oo [ i ]);
1072  }
1073 
1074  * this = tmp;
1075  }
1076 }
1077 
1078 
1079 std :: ostream &operator << ( std :: ostream & out, const FloatArray & x )
1080 {
1081  out << x.giveSize();
1082  for ( double xi : x ) {
1083  out << " " << xi;
1084  }
1085  return out;
1086 }
1087 
1088 //void FloatArray :: beReducedVectorFormOfStrain(const FloatMatrix &aMatrix)
1089 //{
1090 // // Revrites the matrix on vector form (symmetrized matrix used), order: 11, 22, 33, 23, 13, 12
1091 //# ifdef DEBUG
1092 // if ( aMatrix.giveNumberOfColumns() !=3 || aMatrix.giveNumberOfColumns() !=3) {
1093 // OOFEM_ERROR("matrix dimension is not 3x3");
1094 // }
1095 //
1096 //# endif
1097 //
1098 // this->values.resize(6);
1099 // this->at(1) = aMatrix.at(1,1); this->at(2) = aMatrix.at(2,2); this->at(3) = aMatrix.at(3,3);
1100 // this->at(4) = ( aMatrix.at(2,3) + aMatrix.at(3,2) ); // Shear strains multiplied with a factor of 2
1101 // this->at(5) = ( aMatrix.at(1,3) + aMatrix.at(3,1) );
1102 // this->at(6) = ( aMatrix.at(1,2) + aMatrix.at(2,1) );
1103 //}
1104 
1105 void FloatArray :: power(const double exponent)
1106 {
1107  for ( double &x : *this ) {
1108  x = pow(x, exponent);
1109  }
1110 }
1111 
1112 
1113 
1114 void FloatArray :: beColumnOf(const FloatMatrix &mat, int col)
1115 {
1117  mat.copyColumn(*this, col);
1118 }
1119 } // end namespace oofem
int giveNumberOfColumns() const
Returns number of columns of receiver.
Definition: floatmatrix.h:158
void subtract(const FloatArray &src)
Subtracts array src to receiver.
Definition: floatarray.C:258
void checkBounds(int i) const
Checks size of receiver towards requested bounds.
Definition: floatarray.C:132
double & operator()(int i)
Coefficient access function.
Definition: floatarray.h:153
void copySubVector(const FloatArray &src, int si)
Copy the given vector as sub-vector to receiver.
Definition: floatarray.C:864
void beVectorProductOf(const FloatArray &v1, const FloatArray &v2)
Computes vector product (or cross product) of vectors given as parameters, , and stores the result in...
Definition: floatarray.C:415
friend std::ostream & operator<<(std::ostream &out, const FloatArray &x)
Definition: floatarray.C:1079
std::vector< double > values
Stored values.
Definition: floatarray.h:86
FloatArray & operator*=(FloatArray &x, const double &a)
Vector multiplication by scalar.
Definition: floatarray.C:934
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
contextIOResultType storeYourself(DataStream &stream) const
Definition: floatarray.C:872
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
int max(int i, int j)
Returns bigger value form two given decimals.
Definition: mathfem.h:71
int minimum() const
Finds the minimum component in the array.
Definition: intarray.C:184
FloatArray & operator-=(FloatArray &x, const FloatArray &y)
Definition: floatarray.C:974
void beVectorForm(const FloatMatrix &aMatrix)
Reciever will be a vector with 9 components formed from a 3x3 matrix.
Definition: floatarray.C:992
#define FAST_RESIZE(newsize)
Definition: floatarray.C:54
General IO error.
const double * givePointer() const
Exposes the internal values of the matrix.
Definition: floatmatrix.h:558
double product() const
Computes the product of receiver values.
Definition: floatarray.C:858
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
virtual int read(int *data, int count)=0
Reads count integer values into array pointed by data.
void beMaxOf(const FloatArray &a, const FloatArray &b)
Sets receiver to maximum of a or b&#39;s respective elements.
Definition: floatarray.C:288
void changeComponentOrder()
Swaps the fourth and sixth index in the array.
Definition: floatarray.C:1055
void rotatedWith(FloatMatrix &r, char mode)
Returns the receiver a rotated according the change-of-base matrix r.
Definition: floatarray.C:799
bool isFinite() const
Returns true if no element is NAN or infinite.
Definition: floatarray.C:72
void beDifferenceOf(const FloatArray &a, const FloatArray &b)
Sets receiver to be a - b.
Definition: floatarray.C:341
void beColumnOf(const FloatMatrix &mat, int col)
Reciever will be set to a given column in a matrix.
Definition: floatarray.C:1114
void checkSizeTowards(const IntArray &loc)
Checks size of receiver towards values stored in loc array.
Definition: floatarray.C:593
double distance(const FloatArray &x) const
Computes the distance between position represented by receiver and position given as parameter...
Definition: floatarray.C:489
double sum() const
Computes the sum of receiver values.
Definition: floatarray.C:852
void beScaled(double s, const FloatArray &b)
Sets receiver to be .
Definition: floatarray.C:146
void beMinOf(const FloatArray &a, const FloatArray &b)
Sets receiver to be minimum of a or b&#39;s respective elements.
Definition: floatarray.C:315
double dotProduct(const FloatArray &x) const
Computes the dot product (or inner product) of receiver and argument.
Definition: floatarray.C:463
#define OOFEM_FATAL(...)
Macros for printing errors.
Definition: error.h:60
virtual void printYourselfToFile(const std::string filename, const bool showDimensions=true) const
Print receiver to file.
Definition: floatarray.C:770
int givePackSize(DataStream &buff) const
Definition: floatarray.C:921
#define OOFEM_ERROR(...)
Definition: error.h:61
double computeSquaredNorm() const
Computes the square of the norm.
Definition: floatarray.C:846
void addSubVector(const FloatArray &src, int si)
Adds the given vector as sub-vector to receiver.
Definition: floatarray.C:399
contextIOResultType restoreYourself(DataStream &stream)
Definition: floatarray.C:895
int maximum() const
Finds the maximum component in the array.
Definition: intarray.C:195
void power(const double exponent)
Raise each element to its power.
Definition: floatarray.C:1105
FloatArray & operator=(const FloatArray &src)
Assignment operator.
Definition: floatarray.h:111
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
Definition: floatarray.C:676
void resizeWithValues(int s, int allocChunk=0)
Checks size of receiver towards requested bounds.
Definition: floatarray.C:615
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix^T and anArray.
Definition: floatarray.C:708
bool isEmpty() const
Returns true if receiver is empty.
Definition: floatarray.h:222
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
double distance_square(const FloatArray &iP1, const FloatArray &iP2, double &oXi, double &oXiUnbounded) const
Definition: floatarray.C:499
void hardResize(int s)
Resizes the size of the receiver to requested bounds.
Definition: floatarray.C:639
void reserve(int s)
Allocates enough size to fit s, and clears the array.
Definition: floatarray.C:608
int giveIndexMaxElem(void)
Returns index (between 1 and Size) of maximum element in the array.
Definition: floatarray.C:447
int giveIndexMinElem(void)
Returns index (between 1 and Size) of minimum element in the array.
Definition: floatarray.C:431
Class representing vector of real numbers.
Definition: floatarray.h:82
double & operator[](int i)
Definition: floatarray.h:156
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
virtual int givePackSizeOfDouble(int count)=0
void beSubArrayOf(const FloatArray &src, const IntArray &indx)
Extract sub vector form src array and stores the result into receiver.
Definition: floatarray.C:379
void assemble(const FloatArray &fe, const IntArray &loc)
Assembles the array fe (typically, the load vector of a finite element) into the receiver, using loc as location array.
Definition: floatarray.C:551
double norm(const FloatArray &x)
Definition: floatarray.C:985
virtual void pY() const
Print receiver on stdout with high accuracy.
Definition: floatarray.C:787
double computeNorm() const
Computes the norm (or length) of the vector.
Definition: floatarray.C:840
virtual void printYourself() const
Print receiver on stdout.
Definition: floatarray.C:747
void append(const FloatArray &a)
Appends array to reciever.
Definition: floatarray.C:664
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
void copyColumn(FloatArray &dest, int c) const
Fetches the values from the specified column.
Definition: floatmatrix.C:662
void assembleSquared(const FloatArray &fe, const IntArray &loc)
Assembles the array fe with each component squared.
Definition: floatarray.C:572
void times(double s)
Multiplies receiver with scalar.
Definition: floatarray.C:818
const double * givePointer() const
Gives the pointer to the raw data, breaking encapsulation.
Definition: floatarray.h:469
std::vector< double >::iterator end()
Definition: floatarray.h:92
FloatArray operator*(const double &a, const FloatArray &x)
Definition: floatarray.C:940
FloatArray operator-(const FloatArray &x, const FloatArray &y)
Definition: floatarray.C:961
void beSymVectorForm(const FloatMatrix &aMatrix)
Reciever will be a vector with 6 components formed from a 3x3 matrix.
Definition: floatarray.C:1035
int min(int i, int j)
Returns smaller value from two given decimals.
Definition: mathfem.h:59
void beSymVectorFormOfStrain(const FloatMatrix &aMatrix)
Definition: floatarray.C:1014
FloatArray & operator+=(FloatArray &x, const FloatArray &y)
Definition: floatarray.C:968
int giveSize() const
Definition: intarray.h:203
int giveSize() const
Returns the size of receiver.
Definition: floatarray.h:218
bool containsOnlyZeroes() const
Returns nonzero if all coefficients of the receiver are 0, else returns zero.
Definition: floatarray.C:646
the oofem namespace is to define a context or scope in which all oofem names are defined.
void negated()
Switches the sign of every coefficient of receiver.
Definition: floatarray.C:739
double normalize()
Normalizes receiver.
Definition: floatarray.C:828
int giveNumberOfRows() const
Returns number of rows of receiver.
Definition: floatmatrix.h:156
FloatArray operator+(const FloatArray &x, const FloatArray &y)
Definition: floatarray.C:954
void add(const FloatArray &src)
Adds array src to receiver.
Definition: floatarray.C:156
double dot(const FloatArray &x, const FloatArray &y)
Definition: floatarray.C:980
std::vector< double >::iterator begin()
Definition: floatarray.h:91
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631
void plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
Adds the product .
Definition: floatarray.C:226
virtual int givePackSizeOfInt(int count)=0

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:28 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011