]> https://gitweb.dealii.org/ - dealii.git/commitdiff
ARKode with L:d:(Block)Vector 11660/head
authorSebastian Proell <proell@lnm.mw.tum.de>
Mon, 1 Feb 2021 12:13:41 +0000 (13:13 +0100)
committerSebastian Proell <proell@lnm.mw.tum.de>
Mon, 1 Feb 2021 21:40:12 +0000 (22:40 +0100)
doc/news/changes/minor/20210131Proell [new file with mode: 0644]
source/sundials/arkode.cc
tests/sundials/arkode_03.cc

diff --git a/doc/news/changes/minor/20210131Proell b/doc/news/changes/minor/20210131Proell
new file mode 100644 (file)
index 0000000..9a77648
--- /dev/null
@@ -0,0 +1,5 @@
+New: SUNDIALS N_Vector module now directly operates on different deal.II vectors without internally
+creating copies.
+In particular, ARKode can now also be used with LinearAlgebra::distributed::(Block)Vector.
+<br>
+(Sebastian Proell, 2021/01/31)
index 7cf6d94a4d39c8cc09d7e52e3d168b14f5a7ca8a..3d0187c6d55bf071c633dc379f015f7823ae5aa6 100644 (file)
@@ -24,6 +24,8 @@
 #  include <deal.II/base/utilities.h>
 
 #  include <deal.II/lac/block_vector.h>
+#  include <deal.II/lac/la_parallel_block_vector.h>
+#  include <deal.II/lac/la_parallel_vector.h>
 #  ifdef DEAL_II_WITH_TRILINOS
 #    include <deal.II/lac/trilinos_parallel_block_vector.h>
 #    include <deal.II/lac/trilinos_vector.h>
@@ -1016,12 +1018,22 @@ namespace SUNDIALS
   template class ARKode<Vector<double>>;
   template class ARKode<BlockVector<double>>;
 
+  template class ARKode<LinearAlgebra::distributed::Vector<double>>;
+  template class ARKode<LinearAlgebra::distributed::BlockVector<double>>;
+
 #  if DEAL_II_SUNDIALS_VERSION_GTE(4, 0, 0)
   template struct SundialsOperator<Vector<double>>;
   template struct SundialsOperator<BlockVector<double>>;
+  template struct SundialsOperator<LinearAlgebra::distributed::Vector<double>>;
+  template struct SundialsOperator<
+    LinearAlgebra::distributed::BlockVector<double>>;
 
   template struct SundialsPreconditioner<Vector<double>>;
   template struct SundialsPreconditioner<BlockVector<double>>;
+  template struct SundialsPreconditioner<
+    LinearAlgebra::distributed::Vector<double>>;
+  template struct SundialsPreconditioner<
+    LinearAlgebra::distributed::BlockVector<double>>;
 #  endif
 
 #  ifdef DEAL_II_WITH_MPI
@@ -1045,11 +1057,11 @@ namespace SUNDIALS
   template class ARKode<PETScWrappers::MPI::BlockVector>;
 
 #        if DEAL_II_SUNDIALS_VERSION_GTE(4, 0, 0)
-  template class SundialsOperator<PETScWrappers::MPI::Vector>;
-  template class SundialsOperator<PETScWrappers::MPI::BlockVector>;
+  template struct SundialsOperator<PETScWrappers::MPI::Vector>;
+  template struct SundialsOperator<PETScWrappers::MPI::BlockVector>;
 
-  template class SundialsPreconditioner<PETScWrappers::MPI::Vector>;
-  template class SundialsPreconditioner<PETScWrappers::MPI::BlockVector>;
+  template struct SundialsPreconditioner<PETScWrappers::MPI::Vector>;
+  template struct SundialsPreconditioner<PETScWrappers::MPI::BlockVector>;
 #        endif
 #      endif // PETSC_USE_COMPLEX
 #    endif   // DEAL_II_WITH_PETSC
index db0a27ce609d05646caa5f190641fbd91782f849..63f591339e4370acddc4087821903bed01d5cc53 100644 (file)
 #include <deal.II/base/parameter_handler.h>
 
 #include <deal.II/lac/full_matrix.h>
-#include <deal.II/lac/vector.h>
+#include <deal.II/lac/la_parallel_vector.h>
 
 #include <deal.II/sundials/arkode.h>
 
 #include "../tests.h"
 
 
-// Test implicit-explicit time stepper, no jacobian.
+// Test implicit-explicit time stepper, no jacobian. Use L:d:V (in serial)
 // Brusselator benchmark
 
 /**
@@ -54,7 +54,7 @@ main(int argc, char **argv)
   Utilities::MPI::MPI_InitFinalize mpi_initialization(
     argc, argv, numbers::invalid_unsigned_int);
 
-  typedef Vector<double> VectorType;
+  using VectorType = LinearAlgebra::distributed::Vector<double>;
 
   ParameterHandler                             prm;
   SUNDIALS::ARKode<VectorType>::AdditionalData data;
@@ -108,7 +108,7 @@ main(int argc, char **argv)
     return 0;
   };
 
-  Vector<double> y(3);
+  VectorType y(3);
   y[0] = u0;
   y[1] = v0;
   y[2] = w0;

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.