<ol>
-<li> Fixed: FETools::interpolation_difference with PETSc.
+<li> Changed: FEValuesExtractors::Scalar, FEValuesExtractors::Vector and
+FEValuesExtractors::SymmetricTensor could not be default constructed, and
+consequently one could not easily put them into arrays (where they would
+be default constructed when changing the size, and later assigned useful
+values). These classes can now be default constructed to invalid
+values, but can of course not be used in any useful way.
+<br>
+(Wolfgang Bangerth, 2013/03/15)
+
+<li> Fixed: FETools::interpolation_difference did not work with PETSc.
+This is now fixed.
<br>
(Timo Heister, 2013/03/01)
//---------------------------------------------------------------------------
// $Id$
//
-// Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010, 2011, 2012 by the deal.II authors
+// Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010, 2011, 2012, 2013 by the deal.II authors
//
// This file is subject to QPL and may not be distributed
// without copyright and license information. Please refer
*/
unsigned int component;
+ /**
+ * Default constructor. Initialize the
+ * object with an invalid component. This leads to
+ * an object that can not be used, but it allows
+ * objects of this kind to be put into arrays that
+ * require a default constructor upon resizing the
+ * array, and then later assigning a suitable
+ * object to each element of the array.
+ */
+ Scalar ();
+
/**
* Constructor. Take the selected
* vector component as argument.
*/
unsigned int first_vector_component;
+ /**
+ * Default constructor. Initialize the
+ * object with an invalid component. This leads to
+ * an object that can not be used, but it allows
+ * objects of this kind to be put into arrays that
+ * require a default constructor upon resizing the
+ * array, and then later assigning a suitable
+ * object to each element of the array.
+ */
+ Vector ();
+
/**
* Constructor. Take the first
* component of the selected vector
*/
unsigned int first_tensor_component;
+ /**
+ * Default constructor. Initialize the
+ * object with an invalid component. This leads to
+ * an object that can not be used, but it allows
+ * objects of this kind to be put into arrays that
+ * require a default constructor upon resizing the
+ * array, and then later assigning a suitable
+ * object to each element of the array.
+ */
+ SymmetricTensor ();
+
/**
* Constructor. Take the first
* component of the selected tensor
namespace FEValuesExtractors
{
+ inline
+ Scalar::Scalar ()
+ :
+ component (numbers::invalid_unsigned_int)
+ {}
+
+
+
inline
Scalar::Scalar (const unsigned int component)
:
+ inline
+ Vector::Vector ()
+ :
+ first_vector_component (numbers::invalid_unsigned_int)
+ {}
+
+
inline
Vector::Vector (const unsigned int first_vector_component)
:
first_vector_component (first_vector_component)
{}
+
+ template <int rank>
+ inline
+ SymmetricTensor<rank>::SymmetricTensor ()
+ :
+ first_tensor_component(numbers::invalid_unsigned_int)
+ {}
+
+
template <int rank>
inline
SymmetricTensor<rank>::SymmetricTensor (const unsigned int first_tensor_component)
--- /dev/null
+//----------------------------------------------------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 2007, 2008, 2013 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.
+//
+//----------------------------------------------------------------------
+
+
+// make sure FEValuesExtractors::Scalar can be default-constructed but that it
+// produces an invalid, unusable object
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/fe/fe_values.h>
+
+#include <fstream>
+
+
+
+template<int dim>
+void test (const Triangulation<dim>& tr,
+ const FiniteElement<dim>& fe)
+{
+ DoFHandler<dim> dof(tr);
+ dof.distribute_dofs(fe);
+
+ const QGauss<dim> quadrature(2);
+ FEValues<dim> fe_values (fe, quadrature,
+ update_values);
+ fe_values.reinit (dof.begin_active());
+
+ FEValuesExtractors::Scalar extr; // invalid object
+ try
+ {
+ fe_values[extr]; // invalid access
+ }
+ catch (const std::exception &exc)
+ {
+ goto ok;
+ }
+
+ Assert (false, ExcMessage ("No exception!?"));
+
+ ok:
+ ;
+}
+
+
+
+template<int dim>
+void test()
+{
+ Triangulation<dim> tr;
+ GridGenerator::hyper_cube(tr);
+
+ FESystem<dim> fe (FE_Q<dim>(1), 1);
+ test(tr, fe);
+}
+
+
+int main()
+{
+ std::ofstream logfile ("fe_values_view_invalid_01/output");
+ deallog << std::setprecision (2);
+
+ deallog.attach(logfile);
+ deallog.depth_console (0);
+ deallog.threshold_double(1.e-7);
+
+ deal_II_exceptions::disable_abort_on_exception();
+
+ test<1>();
+ test<2>();
+ test<3>();
+}
--- /dev/null
+
+DEAL::ExcIndexRange (scalar.component, 0, fe_values_views_cache.scalars.size())
+DEAL::ExcIndexRange (scalar.component, 0, fe_values_views_cache.scalars.size())
+DEAL::ExcIndexRange (scalar.component, 0, fe_values_views_cache.scalars.size())
--- /dev/null
+//----------------------------------------------------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 2007, 2008, 2013 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.
+//
+//----------------------------------------------------------------------
+
+
+// make sure FEValuesExtractors::Vector can be default-constructed but that it
+// produces an invalid, unusable object
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/fe/fe_values.h>
+
+#include <fstream>
+
+
+
+template<int dim>
+void test (const Triangulation<dim>& tr,
+ const FiniteElement<dim>& fe)
+{
+ DoFHandler<dim> dof(tr);
+ dof.distribute_dofs(fe);
+
+ const QGauss<dim> quadrature(2);
+ FEValues<dim> fe_values (fe, quadrature,
+ update_values);
+ fe_values.reinit (dof.begin_active());
+
+ FEValuesExtractors::Vector extr; // invalid object
+ try
+ {
+ fe_values[extr]; // invalid access
+ }
+ catch (const std::exception &exc)
+ {
+ goto ok;
+ }
+
+ Assert (false, ExcMessage ("No exception!?"));
+
+ ok:
+ ;
+}
+
+
+
+template<int dim>
+void test()
+{
+ Triangulation<dim> tr;
+ GridGenerator::hyper_cube(tr);
+
+ FESystem<dim> fe (FE_Q<dim>(1), dim);
+ test(tr, fe);
+}
+
+
+int main()
+{
+ std::ofstream logfile ("fe_values_view_invalid_02/output");
+ deallog << std::setprecision (2);
+
+ deallog.attach(logfile);
+ deallog.depth_console (0);
+ deallog.threshold_double(1.e-7);
+
+ deal_II_exceptions::disable_abort_on_exception();
+
+ test<1>();
+ test<2>();
+ test<3>();
+}
--- /dev/null
+
+DEAL::ExcIndexRange (vector.first_vector_component, 0, fe_values_views_cache.vectors.size())
+DEAL::ExcIndexRange (vector.first_vector_component, 0, fe_values_views_cache.vectors.size())
+DEAL::ExcIndexRange (vector.first_vector_component, 0, fe_values_views_cache.vectors.size())
--- /dev/null
+//----------------------------------------------------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 2007, 2008, 2013 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.
+//
+//----------------------------------------------------------------------
+
+
+// make sure FEValuesExtractors::SymmetricTensor can be default-constructed
+// but that it produces an invalid, unusable object
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/fe/fe_values.h>
+
+#include <fstream>
+
+
+
+template<int dim>
+void test (const Triangulation<dim>& tr,
+ const FiniteElement<dim>& fe)
+{
+ DoFHandler<dim> dof(tr);
+ dof.distribute_dofs(fe);
+
+ const QGauss<dim> quadrature(2);
+ FEValues<dim> fe_values (fe, quadrature,
+ update_values);
+ fe_values.reinit (dof.begin_active());
+
+ FEValuesExtractors::SymmetricTensor<2> extr; // invalid object
+ try
+ {
+ fe_values[extr]; // invalid access
+ }
+ catch (const std::exception &exc)
+ {
+ goto ok;
+ }
+
+ Assert (false, ExcMessage ("No exception!?"));
+
+ ok:
+ ;
+}
+
+
+
+template<int dim>
+void test()
+{
+ Triangulation<dim> tr;
+ GridGenerator::hyper_cube(tr);
+
+ FESystem<dim> fe (FE_Q<dim>(1), dim*(dim+1)/2);
+ test(tr, fe);
+}
+
+
+int main()
+{
+ std::ofstream logfile ("fe_values_view_invalid_03/output");
+ deallog << std::setprecision (2);
+
+ deallog.attach(logfile);
+ deallog.depth_console (0);
+ deallog.threshold_double(1.e-7);
+
+ deal_II_exceptions::disable_abort_on_exception();
+
+ test<1>();
+ test<2>();
+ test<3>();
+}
--- /dev/null
+
+DEAL::ExcIndexRange (tensor.first_tensor_component, 0, fe_values_views_cache.symmetric_second_order_tensors.size())
+DEAL::ExcIndexRange (tensor.first_tensor_component, 0, fe_values_views_cache.symmetric_second_order_tensors.size())
+DEAL::ExcIndexRange (tensor.first_tensor_component, 0, fe_values_views_cache.symmetric_second_order_tensors.size())