- Namespace deal_II_numbers.
- Triangulation::distort_random.
- Triangulation::clear_user_pointers.
+- The deprecated constructor of SparseILU.
<!-- ----------- GENERAL IMPROVEMENTS ----------------- -->
*/
SparseILU ();
- /**
- * @deprecated This method is deprecated, and left for backward
- * compatibility. It will be removed in later versions. Instead, pass the
- * sparsity pattern that you want used for the decomposition in the
- * AdditionalData structure.
- */
- SparseILU (const SparsityPattern &sparsity) DEAL_II_DEPRECATED;
-
/**
* Make SparseLUDecomposition::AdditionalData accessible to this class as
* well.
-template <typename number>
-SparseILU<number>::SparseILU (const SparsityPattern &sparsity)
-{
- SparseMatrix<number>::reinit(sparsity);
-}
-
-
-
-
template <typename number>
template <typename somenumber>
void SparseILU<number>::initialize (const SparseMatrix<somenumber> &matrix,
// ---------------------------------------------------------------------
//
-// Copyright (C) 2007 - 2013 by the deal.II authors
+// Copyright (C) 2007 - 2013, 2015 by the deal.II authors
//
// This file is part of the deal.II library.
//
rhs(1) = rhs(3) = 0.0975;
- SparseILU<double> ilu (sparsity_pattern);
- ilu.decompose (M);
+ SparseILU<double>::AdditionalData data;
+ data.use_this_sparsity = &sparsity_pattern;
+ SparseILU<double> ilu;
+ ilu.initialize (M, data);
Vector<double> solution (4);
// this would fail with elements of
// ---------------------------------------------------------------------
//
-// Copyright (C) 2001 - 2013 by the deal.II authors
+// Copyright (C) 2001 - 2013, 2015 by the deal.II authors
//
// This file is part of the deal.II library.
//
Assert (false, ExcNotImplemented());
};
ilu_pattern.compress();
- SparseILU<double> ilu (ilu_pattern);
- ilu.decompose (A);
+
+ SparseILU<double>::AdditionalData data;
+ data.use_this_sparsity = &ilu_pattern;
+ SparseILU<double> ilu;
+ ilu.initialize (A, data);
// now for three test vectors v
// determine norm of
// ---------------------------------------------------------------------
//
-// Copyright (C) 2001 - 2013 by the deal.II authors
+// Copyright (C) 2001 - 2013, 2015 by the deal.II authors
//
// This file is part of the deal.II library.
//
Assert (false, ExcNotImplemented());
};
ilu_pattern.compress();
- SparseILU<double> ilu (ilu_pattern);
- ilu.decompose (A);
+ SparseILU<double>::AdditionalData data;
+ data.use_this_sparsity = &ilu_pattern;
+ SparseILU<double> ilu;
+ ilu.initialize (A, data);
// now for three test vectors v
// determine norm of