* transformation F. That is to say if $DF$ is square, it computes
* $\det(DF)$, in case DF is not square returns $\sqrt{\det(DF^T * DF)}$.
*/
- double determinant () const;
+ Number determinant () const;
/**
* Assuming that the current object stores the Jacobian of a mapping
template <int order, int dim, int spacedim, typename Number>
inline
-double
-DerivativeForm<order,dim,spacedim,Number>::determinant () const
+Number DerivativeForm<order, dim, spacedim, Number>::determinant() const
{
Assert( order==1, ExcMessage("Only for order == 1."));
if (dim == spacedim)
else
{
Assert( spacedim>dim, ExcMessage("Only for spacedim>dim."));
- const DerivativeForm<1,spacedim,dim> DF_t = this->transpose();
+ const DerivativeForm<1,spacedim,dim,Number> DF_t = this->transpose();
Tensor<2,dim,Number> G; //First fundamental form
for (unsigned int i=0; i<dim; ++i)
for (unsigned int j=0; j<dim; ++j)
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2000 - 2016 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+// exercise DerivativeForm::norm() for complex data types
+
+
+#include "../tests.h"
+#include <deal.II/base/derivative_form.h>
+
+
+template<int dim, int spacedim>
+void test()
+{
+ DerivativeForm<1,dim,spacedim,std::complex<double> > dF;
+ for (unsigned int i=0; i<spacedim; ++i)
+ for (unsigned int j=0; j<dim; ++j)
+ {
+ dF[i][j] = std::complex<double>(i+2*j+1,i+2*j+1);
+ }
+
+ // output the determinants of these objects
+ deallog << "det(dF): " << dF.determinant() << std::endl;
+}
+
+int main()
+{
+ initlog();
+ test<1,1>();
+ test<1,2>();
+ test<1,3>();
+ test<2,2>();
+ test<2,3>();
+ test<3,3>();
+}
--- /dev/null
+
+DEAL::det(dF): (1.00000,1.00000)
+DEAL::det(dF): (2.23607,2.23607)
+DEAL::det(dF): (3.74166,3.74166)
+DEAL::det(dF): (0.00000,-4.00000)
+DEAL::det(dF): (5.99952e-16,9.79796)
+DEAL::det(dF): (0.00000,0.00000)