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.); \ 61 #ifdef __LAPACK_MODULE 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);
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);
75 if(!std::isfinite(val)) {
150 for (
int i = 0; i < this->
giveSize(); ++i ) {
151 (*this) [ i ] = s * b [ i ];
177 #ifdef __LAPACK_MODULE 183 for (
int i = 0; i < this->
giveSize(); i++ ) {
184 (*this) [ i ] += b [ i ];
192 for (
double &x: *
this ) {
214 #ifdef __LAPACK_MODULE 219 for (
int i = 0; i < this->
giveSize(); ++i ) {
220 (*this) [ i ] += factor * b [ i ];
233 this->
values.assign( nColumns, 0. );
242 #ifdef __LAPACK_MODULE 245 dgemv_(
"t", & nRows, & nColumns, & dV, b.
givePointer(), & nRows, s.
givePointer(), & inc, & beta, this->
givePointer(), & inc, nColumns, nColumns, nRows);
247 for (
int i = 1; i <= nColumns; i++ ) {
249 for (
int j = 1; j <= nRows; j++ ) {
250 sum += b.
at(j, i) * s.
at(j);
252 this->
at(i) += sum * dV;
268 for (
int i = 0; i < this->
giveSize(); ++i ) {
269 (*this) [ i ] = -src [ i ];
282 for (
int i = 0; i < this->
giveSize(); ++i ) {
283 (*this) [ i ] -= src [ i ];
309 for (
int i = 0; i < n; i++ ) {
310 (*this) [ i ] =
max( a [ i ], b [ i ] );
335 for (
int i = 0; i < n; i++ ) {
336 (*this) [ i ] =
min( a [ i ], b [ i ] );
351 for (
int i = 0; i < this->
giveSize(); ++i ) {
352 (*this) [ i ] = a [ i ] - b [ i ];
357 for (
int i = 0; i < a.
giveSize(); ++i ) {
358 this->
values.push_back( a[i] - b[i] );
373 for (
int i = 0; i < n; ++i ) {
374 (*this) [ i ] = a [ i ] - b [ i ];
393 for (
int i = 1; i <= n; i++ ) {
394 this->
at(i) = src.
at( indx.
at(i) );
409 for (
int i = 0; i < n; i++ ) {
410 (*this) [si + i] += src [ i ];
420 OOFEM_ERROR(
"size mismatch, size is not equal to 3");
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);
437 double val = (*this) [ 0 ];
438 for (
int i = 1; i < this->
giveSize(); i++ ) {
439 if ( val > (*
this) [ i ] ) {
453 double val = (*this) [ 0 ];
454 for (
int i = 1; i < this->
giveSize(); i++ ) {
455 if ( val < (*
this) [ i ] ) {
472 return std::inner_product(this->
begin(), this->
end(), x.
begin(), 0.);
485 return std::inner_product(this->
begin(), this->
begin()+size, x.
begin(), 0.);
505 const double s = (
dot(*
this, iP2) -
dot(*
this, iP1) ) - (
dot(iP1, iP2) -
dot(iP1, iP1) );
521 const FloatArray q = ( 1.0 - oXi ) * iP1 + oXi * iP2;
542 for (
int i = 1; i <= s; ++i ) {
543 double dx = this->
at(i) - from.at(i);
563 for (
int i = 1; i <= n; i++ ) {
566 this->
at(ii) += fe.
at(i);
584 for (
int i = 1; i <= n; i++ ) {
587 this->
at(ii) += fe.
at(i) * fe.
at(i);
598 for (
int p : loc ) {
599 high =
max( high, p );
603 this->
values.resize(high);
618 if ( allocChunk < 0 ) {
619 OOFEM_FATAL(
"allocChunk must be non-negative; %d", allocChunk);
624 if ( allocChunk > 0 && (
int)this->
values.capacity() < n ) {
625 this->
values.reserve(n + allocChunk);
641 this->
values.assign(n, 0.);
642 this->
values.shrink_to_fit();
648 for (
double x : *
this ) {
660 std::fill(this->
begin(), this->
end(), 0.);
672 this->
values.push_back(a);
691 #ifdef __LAPACK_MODULE 692 double alpha = 1., beta = 0.;
694 dgemv_(
"n", & nRows, & nColumns, & alpha, aMatrix.
givePointer(), & nRows, anArray.
givePointer(), & inc, & beta, this->
givePointer(), & inc, nColumns, nColumns, nRows);
696 for (
int i = 1; i <= nRows; i++ ) {
698 for (
int j = 1; j <= nColumns; j++ ) {
699 sum += aMatrix.
at(i, j) * anArray.
at(j);
722 #ifdef __LAPACK_MODULE 723 double alpha = 1., beta = 0.;
725 dgemv_(
"t", & nRows, & nColumns, & alpha, aMatrix.
givePointer(), & nRows, anArray.
givePointer(), & inc, & beta, this->
givePointer(), & inc, nColumns, nColumns, nRows);
727 for (
int i = 1; i <= nColumns; i++ ) {
729 for (
int j = 1; j <= nRows; j++ ) {
730 sum += aMatrix.
at(j, i) * anArray.
at(j);
741 for (
double &x: *
this ) {
750 printf(
"FloatArray of size : %d \n", this->
giveSize());
751 for (
double x: *
this ) {
752 printf(
"%10.3e ", x );
762 printf(
"%s (%d): \n", name.c_str(), this->
giveSize());
763 for (
double x: *
this ) {
764 printf(
"%10.3e ", x );
773 std :: ofstream arrayfile (filename);
774 if (arrayfile.is_open()) {
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";
791 for (
double x: *
this ) {
792 printf(
"%20.14e; ", x );
808 }
else if ( mode ==
'n' ) {
822 for (
double &x : *
this ) {
831 if ( norm < 1.e-80 ) {
832 OOFEM_ERROR(
"cannot norm receiver, norm is too small");
835 this->
times(1. / norm);
848 return std::inner_product(this->
begin(), this->
end(), this->
begin(), 0.);
854 return std::accumulate(this->
begin(), this->
end(), 0.);
860 return std::accumulate(this->
begin(), this->
end(), 1.0, [](
double a,
double b) {
return a*b; });
880 if ( !stream.write(size) ) {
886 if ( !stream.write(this->givePointer(), size) ) {
903 if ( !stream.
read(size) ) {
907 this->
values.resize(size);
911 if ( !stream.
read(this->givePointer(), size) ) {
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) )
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) )
1061 std :: swap( this->
at(4), this->
at(6) );
1062 }
else if ( this->
giveSize() == 9 ) {
1065 const int abq2oo [ 9 ] = {
1066 1, 2, 3, 6, 5, 4, 7, 9, 8
1070 for (
int i = 0; i < 9; i++ ) {
1071 tmp(i) = this->
at(abq2oo [ i ]);
1082 for (
double xi : x ) {
1107 for (
double &x : *
this ) {
1108 x = pow(x, exponent);
int giveNumberOfColumns() const
Returns number of columns of receiver.
void subtract(const FloatArray &src)
Subtracts array src to receiver.
void checkBounds(int i) const
Checks size of receiver towards requested bounds.
double & operator()(int i)
Coefficient access function.
void copySubVector(const FloatArray &src, int si)
Copy the given vector as sub-vector to receiver.
void beVectorProductOf(const FloatArray &v1, const FloatArray &v2)
Computes vector product (or cross product) of vectors given as parameters, , and stores the result in...
friend std::ostream & operator<<(std::ostream &out, const FloatArray &x)
std::vector< double > values
Stored values.
FloatArray & operator*=(FloatArray &x, const double &a)
Vector multiplication by scalar.
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
contextIOResultType storeYourself(DataStream &stream) const
double & at(int i)
Coefficient access function.
int max(int i, int j)
Returns bigger value form two given decimals.
int minimum() const
Finds the minimum component in the array.
FloatArray & operator-=(FloatArray &x, const FloatArray &y)
void beVectorForm(const FloatMatrix &aMatrix)
Reciever will be a vector with 9 components formed from a 3x3 matrix.
#define FAST_RESIZE(newsize)
const double * givePointer() const
Exposes the internal values of the matrix.
double product() const
Computes the product of receiver values.
Class implementing an array of integers.
int & at(int i)
Coefficient access function.
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's respective elements.
void changeComponentOrder()
Swaps the fourth and sixth index in the array.
void rotatedWith(FloatMatrix &r, char mode)
Returns the receiver a rotated according the change-of-base matrix r.
bool isFinite() const
Returns true if no element is NAN or infinite.
void beDifferenceOf(const FloatArray &a, const FloatArray &b)
Sets receiver to be a - b.
void beColumnOf(const FloatMatrix &mat, int col)
Reciever will be set to a given column in a matrix.
void checkSizeTowards(const IntArray &loc)
Checks size of receiver towards values stored in loc array.
double distance(const FloatArray &x) const
Computes the distance between position represented by receiver and position given as parameter...
double sum() const
Computes the sum of receiver values.
void beScaled(double s, const FloatArray &b)
Sets receiver to be .
void beMinOf(const FloatArray &a, const FloatArray &b)
Sets receiver to be minimum of a or b's respective elements.
double dotProduct(const FloatArray &x) const
Computes the dot product (or inner product) of receiver and argument.
#define OOFEM_FATAL(...)
Macros for printing errors.
virtual void printYourselfToFile(const std::string filename, const bool showDimensions=true) const
Print receiver to file.
int givePackSize(DataStream &buff) const
double computeSquaredNorm() const
Computes the square of the norm.
void addSubVector(const FloatArray &src, int si)
Adds the given vector as sub-vector to receiver.
contextIOResultType restoreYourself(DataStream &stream)
int maximum() const
Finds the maximum component in the array.
void power(const double exponent)
Raise each element to its power.
FloatArray & operator=(const FloatArray &src)
Assignment operator.
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
void resizeWithValues(int s, int allocChunk=0)
Checks size of receiver towards requested bounds.
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix^T and anArray.
bool isEmpty() const
Returns true if receiver is empty.
double at(int i, int j) const
Coefficient access function.
double distance_square(const FloatArray &iP1, const FloatArray &iP2, double &oXi, double &oXiUnbounded) const
void hardResize(int s)
Resizes the size of the receiver to requested bounds.
void reserve(int s)
Allocates enough size to fit s, and clears the array.
int giveIndexMaxElem(void)
Returns index (between 1 and Size) of maximum element in the array.
int giveIndexMinElem(void)
Returns index (between 1 and Size) of minimum element in the array.
Class representing vector of real numbers.
double & operator[](int i)
Implementation of matrix containing floating point numbers.
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.
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.
double norm(const FloatArray &x)
virtual void pY() const
Print receiver on stdout with high accuracy.
double computeNorm() const
Computes the norm (or length) of the vector.
virtual void printYourself() const
Print receiver on stdout.
void append(const FloatArray &a)
Appends array to reciever.
void zero()
Zeroes all coefficients of receiver.
void copyColumn(FloatArray &dest, int c) const
Fetches the values from the specified column.
void assembleSquared(const FloatArray &fe, const IntArray &loc)
Assembles the array fe with each component squared.
void times(double s)
Multiplies receiver with scalar.
const double * givePointer() const
Gives the pointer to the raw data, breaking encapsulation.
std::vector< double >::iterator end()
FloatArray operator*(const double &a, const FloatArray &x)
FloatArray operator-(const FloatArray &x, const FloatArray &y)
void beSymVectorForm(const FloatMatrix &aMatrix)
Reciever will be a vector with 6 components formed from a 3x3 matrix.
int min(int i, int j)
Returns smaller value from two given decimals.
void beSymVectorFormOfStrain(const FloatMatrix &aMatrix)
FloatArray & operator+=(FloatArray &x, const FloatArray &y)
int giveSize() const
Returns the size of receiver.
bool containsOnlyZeroes() const
Returns nonzero if all coefficients of the receiver are 0, else returns zero.
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.
double normalize()
Normalizes receiver.
int giveNumberOfRows() const
Returns number of rows of receiver.
FloatArray operator+(const FloatArray &x, const FloatArray &y)
void add(const FloatArray &src)
Adds array src to receiver.
double dot(const FloatArray &x, const FloatArray &y)
std::vector< double >::iterator begin()
void resize(int s)
Resizes receiver towards requested size.
void plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
Adds the product .
virtual int givePackSizeOfInt(int count)=0