// $Id$
// Version: $Name$
//
-// Copyright (C) 2005, 2006 by the deal.II authors
+// Copyright (C) 2005, 2006, 2008 by the deal.II authors
//
// This file is subject to QPL and may not be distributed
// without copyright and license information. Please refer
static const unsigned int
n_rank2_components = (dim*dim + dim)/2;
+ /**
+ * Number of independent components of a
+ * symmetric tensor of rank 4.
+ */
+ static const unsigned int
+ n_independent_components = (n_rank2_components *
+ StorageType<2,dim>::n_independent_components);
+
/**
* Declare the type in which we
* actually store the data. Symmetric
* data types.
*/
static const unsigned int dimension = dim;
+
+ /**
+ * An integer denoting the number of
+ * independent components that fully
+ * describe a symmetric tensor. In $d$
+ * space dimensions, this number equals
+ * $\frac 12 (d^2+d)$ for symmetric
+ * tensors of rank 2.
+ */
+ static const unsigned int n_independent_components
+ = internal::SymmetricTensorAccessors::StorageType<rank,dim>::
+ n_independent_components;
/**
* Default constructor. Creates a zero
*/
SymmetricTensor (const Tensor<2,dim> &t);
+ /**
+ * A constructor that creates a symmetric
+ * tensor from an array holding its
+ * independent elements. Using this
+ * constructor assumes that the caller
+ * knows the order in which elements are
+ * stored in symmetric tensors; its use
+ * is therefore discouraged.
+ *
+ * This constructor is currently inly
+ * implemented for symmetric tensors of
+ * rank 2.
+ */
+ SymmetricTensor (const double (&array) [n_independent_components]);
+
/**
* Assignment operator.
*/
}
+template <int rank, int dim>
+inline
+SymmetricTensor<rank,dim>::SymmetricTensor (const double (&array) [n_independent_components])
+ :
+ data (array)
+{}
+
+
template <int rank, int dim>
inline
SymmetricTensor<rank,dim> &
}
+
+/**
+ * Output operator for symmetric tensors of rank 2. Print the elements
+ * consecutively, with a space in between, two spaces between rank 1
+ * subtensors, three between rank 2 and so on. No special amends are made to
+ * represents the symmetry in the output, for example by outputting only the
+ * unique entries.
+ *
+ * @relates SymmetricTensor
+ */
+template <int dim>
+inline
+std::ostream & operator << (std::ostream &out,
+ const SymmetricTensor<2,dim> &t)
+{
+ //make out lives a bit simpler by outputing
+ //the tensor through the operator for the
+ //general Tensor class
+ Tensor<2,dim> tt;
+
+ for (unsigned int i=0; i<dim; ++i)
+ for (unsigned int j=0; j<dim; ++j)
+ tt[i][j] = t[i][j];
+
+ return out << tt;
+}
+
+
+
+/**
+ * Output operator for symmetric tensors of rank 4. Print the elements
+ * consecutively, with a space in between, two spaces between rank 1
+ * subtensors, three between rank 2 and so on. No special amends are made to
+ * represents the symmetry in the output, for example by outputting only the
+ * unique entries.
+ *
+ * @relates SymmetricTensor
+ */
+template <int dim>
+inline
+std::ostream & operator << (std::ostream &out,
+ const SymmetricTensor<4,dim> &t)
+{
+ //make out lives a bit simpler by outputing
+ //the tensor through the operator for the
+ //general Tensor class
+ Tensor<4,dim> tt;
+
+ for (unsigned int i=0; i<dim; ++i)
+ for (unsigned int j=0; j<dim; ++j)
+ for (unsigned int k=0; k<dim; ++k)
+ for (unsigned int l=0; l<dim; ++l)
+ tt[i][j][k][l] = t[i][j][k][l];
+
+ return out << tt;
+}
+
+
DEAL_II_NAMESPACE_CLOSE
#endif
<ol>
+<li> <p> New: The SymmetricTensor class now has a constructor that creates
+an object from an array consisting its independent components.
+<br>
+(WB 2008/02/21)
+</p> </li>
+
+<li> <p> New: There are now output operators (i.e. <code>operator@<@<</code>)
+for the SymmetricTensor class.
+<br>
+(WB 2008/02/12)
+</p> </li>
+
<li> <p> Improved: LogStream is now thread safe and output lines from different
threads are separated now. Additionally, a function LogStream::log_thread_id()
has been added to log the id of the thread printing the message. The semantics
// $Id$
// Version: $Name$
//
-// Copyright (C) 2006 by the deal.II authors
+// Copyright (C) 2006, 2008 by the deal.II authors
//
// This file is subject to QPL and may not be distributed
// without copyright and license information. Please refer
//
//---------------------------- symmetric_tensor_21.cc ---------------------------
-// check SymmetricTensor<2,dim>::operator= (double)
+// check SymmetricTensor<4,dim>::operator= (double)
#include "../tests.h"
#include <base/symmetric_tensor.h>
--- /dev/null
+//---------------------------- symmetric_tensor_22.cc ---------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 2006, 2008 by the deal.II authors
+//
+// This file is subject to QPL and may not be distributed
+// without copyright and license information. Please refer
+// to the file deal.II/doc/license.html for the text and
+// further information on this license.
+//
+//---------------------------- symmetric_tensor_22.cc ---------------------------
+
+// check operator<< for SymmetricTensor<2,dim> and SymmetricTensor<4,dim>
+
+#include "../tests.h"
+#include <base/symmetric_tensor.h>
+#include <base/logstream.h>
+#include <lac/vector.h>
+#include <fstream>
+#include <iomanip>
+
+int main ()
+{
+ std::ofstream logfile("symmetric_tensor_22/output");
+ deallog << std::setprecision(3);
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+ deallog.threshold_double(1.e-10);
+
+ {
+ const unsigned int dim=2;
+ SymmetricTensor<2,dim> t;
+ t[0][0] = 1;
+ t[1][1] = 2;
+ t[0][1] = 3;
+
+ deallog << t << std::endl;
+ }
+
+ {
+ const unsigned int dim=3;
+ SymmetricTensor<4,dim> t;
+ t[0][0][0][0] = t[1][0][1][0] = t[1][1][1][1]
+ = t[2][2][2][2] = t[2][0][2][0] = 3;
+
+ deallog << t << std::endl;
+ }
+
+ deallog << "OK" << std::endl;
+}
--- /dev/null
+
+DEAL::1.00 3.00 3.00 2.00
+DEAL::3.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 3.00 0.00 3.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 3.00 0.00 0.00 0.00 3.00 0.00 0.00 0.00 3.00 0.00 3.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 3.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 3.00 0.00 0.00 0.00 3.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 3.00
+DEAL::OK
--- /dev/null
+//---------------------------- symmetric_tensor_23.cc ---------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 2006, 2008 by the deal.II authors
+//
+// This file is subject to QPL and may not be distributed
+// without copyright and license information. Please refer
+// to the file deal.II/doc/license.html for the text and
+// further information on this license.
+//
+//---------------------------- symmetric_tensor_23.cc ---------------------------
+
+// check operator<< for SymmetricTensor<2,dim> and SymmetricTensor<4,dim>
+
+#include "../tests.h"
+#include <base/symmetric_tensor.h>
+#include <base/logstream.h>
+#include <lac/vector.h>
+#include <fstream>
+#include <iomanip>
+
+int main ()
+{
+ std::ofstream logfile("symmetric_tensor_23/output");
+ deallog << std::setprecision(3);
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+ deallog.threshold_double(1.e-10);
+
+ {
+ SymmetricTensor<2,1> t;
+ t[0][0] = 1;
+
+ double x[1] = { 1 };
+ Assert ((t == SymmetricTensor<2,1>(x)),
+ ExcInternalError());
+ }
+
+ {
+ SymmetricTensor<2,2> t;
+ t[0][0] = 1;
+ t[1][1] = 2;
+ t[0][1] = 3;
+
+ double x[3] = { 1, 2, 3 };
+ Assert ((t == SymmetricTensor<2,2>(x)),
+ ExcInternalError());
+ }
+
+ {
+ SymmetricTensor<2,3> t;
+ t[0][0] = 1;
+ t[1][1] = 2;
+ t[2][2] = 3;
+ t[0][1] = 4;
+ t[0][2] = 5;
+ t[1][2] = 6;
+
+ double x[6] = { 1, 2, 3, 4, 5, 6 };
+ Assert ((t == SymmetricTensor<2,3>(x)),
+ ExcInternalError());
+ }
+
+ deallog << "OK" << std::endl;
+}
--- /dev/null
+
+DEAL::OK