]> https://gitweb.dealii.org/ - dealii.git/commitdiff
extend do_project() for complex numbers with ExcNotImplemented() for solver
authorDenis Davydov <davydden@gmail.com>
Sat, 27 Feb 2016 20:12:59 +0000 (21:12 +0100)
committerDenis Davydov <davydden@gmail.com>
Sat, 27 Feb 2016 21:29:09 +0000 (22:29 +0100)
include/deal.II/numerics/vector_tools.templates.h

index af58c55d2259a1728551244bdc75c0e86f95842f..5f62a6841c5b96a6dc0930c65b3bc5924a92b957 100644 (file)
@@ -752,6 +752,32 @@ namespace VectorTools
       return true;
     }
 
+    template <typename number>
+    void invert_mass_matrix(const SparseMatrix<number> &mass_matrix,
+                            const Vector<number> &tmp,
+                            Vector<number> &vec)
+    {
+      // Allow for a maximum of 5*n steps to reduce the residual by 10^-12. n
+      // steps may not be sufficient, since roundoff errors may accumulate for
+      // badly conditioned matrices
+      ReductionControl      control(5*tmp.size(), 0., 1e-12, false, false);
+      GrowingVectorMemory<Vector<number> > memory;
+      SolverCG<Vector<number> >            cg(control,memory);
+
+      PreconditionSSOR<SparseMatrix<number> > prec;
+      prec.initialize(mass_matrix, 1.2);
+
+      cg.solve (mass_matrix, vec, tmp, prec);
+    }
+
+    template <typename number>
+    void invert_mass_matrix(const SparseMatrix<std::complex<number> > &mass_matrix,
+                            const Vector<std::complex<number> > &tmp,
+                            Vector<std::complex<number> > &vec)
+    {
+      Assert(false, ExcNotImplemented());
+    }
+
 
     /**
      * Generic implementation of the project() function
@@ -825,18 +851,7 @@ namespace VectorTools
           constraints.condense(mass_matrix, tmp);
         }
 
-      // Allow for a maximum of 5*n steps to reduce the residual by 10^-12. n
-      // steps may not be sufficient, since roundoff errors may accumulate for
-      // badly conditioned matrices
-      ReductionControl      control(5*tmp.size(), 0., 1e-12, false, false);
-      GrowingVectorMemory<Vector<number> > memory;
-      SolverCG<Vector<number> >            cg(control,memory);
-
-      PreconditionSSOR<SparseMatrix<number> > prec;
-      prec.initialize(mass_matrix, 1.2);
-
-      cg.solve (mass_matrix, vec, tmp, prec);
-
+      invert_mass_matrix(mass_matrix,tmp,vec);
       constraints.distribute (vec);
 
       // copy vec into vec_result. we can't use vec_result itself above, since

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.