Skip to content
Snippets Groups Projects
Commit 4f69c98c authored by Oliver Sander's avatar Oliver Sander Committed by sander@FU-BERLIN.DE
Browse files

add methods axpy and left-multiplication by a FieldMatrix

[[Imported from SVN: r7135]]
parent f053f105
No related branches found
No related tags found
No related merge requests found
...@@ -35,6 +35,13 @@ class Tensor3 ...@@ -35,6 +35,13 @@ class Tensor3
return norm; return norm;
} }
Tensor3<T,N1,N2,N3>& axpy(const T& alpha, const Tensor3<T,N1,N2,N3>& other)
{
for (int i=0; i<N1; i++)
(*this)[i].axpy(alpha,other[i]);
return *this;
}
static Tensor3<T,N1,N2,N3> product(const Dune::FieldVector<T,N1>& a, const Dune::FieldVector<T,N2>& b, const Dune::FieldVector<T,N3>& c) static Tensor3<T,N1,N2,N3> product(const Dune::FieldVector<T,N1>& a, const Dune::FieldVector<T,N2>& b, const Dune::FieldVector<T,N3>& c)
{ {
Tensor3<T,N1,N2,N3> result; Tensor3<T,N1,N2,N3> result;
...@@ -87,6 +94,22 @@ class Tensor3 ...@@ -87,6 +94,22 @@ class Tensor3
return result; return result;
} }
template <int N4>
friend Tensor3<T,N1,N3,N4> operator*(const Dune::FieldMatrix<T,N1,N2>& a, const Tensor3<T,N2,N3,N4>& b)
{
Tensor3<T,N1,N3,N4> result;
for (int i=0; i<N1; i++)
for (int j=0; j<N3; j++)
for (int k=0; k<N4; k++) {
result[i][j][k] = 0;
for (int l=0; l<N2; l++)
result[i][j][k] += a[i][l]*b[l][j][k];
}
return result;
}
friend Tensor3<T,N1,N2,N3> operator+(const Tensor3<T,N1,N2,N3>& a, const Tensor3<T,N1,N2,N3>& b) friend Tensor3<T,N1,N2,N3> operator+(const Tensor3<T,N1,N2,N3>& a, const Tensor3<T,N1,N2,N3>& b)
{ {
Tensor3<T,N1,N2,N3> result; Tensor3<T,N1,N2,N3> result;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment