Browse Source

Pardiso option

master
Trevor Irons 6 years ago
parent
commit
ac623cc71f
4 changed files with 22 additions and 15 deletions
  1. 7
    3
      CMakeLists.txt
  2. 3
    3
      examples/EMSchur3D-vtk.cpp
  3. 3
    3
      examples/EMSchur3D.cpp
  4. 9
    6
      include/EMSchur3D.h

+ 7
- 3
CMakeLists.txt View File

48
     	CXX_STANDARD_REQUIRED_ON
48
     	CXX_STANDARD_REQUIRED_ON
49
 	)
49
 	)
50
 
50
 
51
-        # No OpenMP
51
+	option ( EMSCHUR3D_ON_KINGSPEAK FALSE )
52
+	if ( EMSCHUR3D_ON_KINGSPEAK )
53
+    # No OpenMP
54
+		add_compile_options(-DHAVE_PARDISO)
52
 	#set (MKLLINK  "-L$MKLROOT/lib/intel64 -lmkl_intel_lp64 -lmkl_sequential -lmkl_core -lpthread -Wl,-rpath=$MKLROOT/lib/intel64")
55
 	#set (MKLLINK  "-L$MKLROOT/lib/intel64 -lmkl_intel_lp64 -lmkl_sequential -lmkl_core -lpthread -Wl,-rpath=$MKLROOT/lib/intel64")
53
-        # OpenMP 
56
+        # OpenMP R3
54
         set (MKLLINK  "-L$MKLROOT/lib/intel64 -Wl,-rpath=$MKLROOT/lib/intel64 -lmkl_intel_lp64 -lmkl_intel_thread -lmkl_core -liomp5 -lpthread")
57
         set (MKLLINK  "-L$MKLROOT/lib/intel64 -Wl,-rpath=$MKLROOT/lib/intel64 -lmkl_intel_lp64 -lmkl_intel_thread -lmkl_core -liomp5 -lpthread")
55
-	target_link_libraries(emschur3d "lemmacore" ${MKLLINK})
58
+		target_link_libraries(emschur3d "lemmacore" ${MKLLINK})
59
+	endif()
56
 
60
 
57
 	# Linking
61
 	# Linking
58
 	target_link_libraries(emschur3d "lemmacore" "fdem1d" )
62
 	target_link_libraries(emschur3d "lemmacore" "fdem1d" )

+ 3
- 3
examples/EMSchur3D-vtk.cpp View File

114
     // And solve
114
     // And solve
115
 
115
 
116
     // Use BiCGSTAB Diagonal preconditioner
116
     // Use BiCGSTAB Diagonal preconditioner
117
-    //auto EM3D = EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::IncompleteLUT<Complex> > >::NewSP();
117
+    auto EM3D = EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::IncompleteLUT<Complex> > >::NewSP();
118
     //auto EM3D = EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::NewSP();
118
     //auto EM3D = EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::NewSP();
119
-    
119
+
120
     //auto EM3D = EMSchur3D< Eigen::PardisoLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::NewSP();
120
     //auto EM3D = EMSchur3D< Eigen::PardisoLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::NewSP();
121
-    auto EM3D = EMSchur3D< Eigen::PardisoLDLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Symmetric > >::NewSP();
121
+    //auto EM3D = EMSchur3D< Eigen::PardisoLDLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Symmetric > >::NewSP();
122
 
122
 
123
     // LS CG
123
     // LS CG
124
     //auto EM3D = EMSchur3D< Eigen::LeastSquaresConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::NewSP();
124
     //auto EM3D = EMSchur3D< Eigen::LeastSquaresConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::NewSP();

+ 3
- 3
examples/EMSchur3D.cpp View File

34
 
34
 
35
     // SUPERLU
35
     // SUPERLU
36
     //auto EM3D = EMSchur3D< Eigen::SuperLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::NewSP();
36
     //auto EM3D = EMSchur3D< Eigen::SuperLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::NewSP();
37
-    
37
+
38
     //auto EM3D = EMSchur3D< Eigen::PardisoLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::NewSP();
38
     //auto EM3D = EMSchur3D< Eigen::PardisoLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::NewSP();
39
-    auto EM3D = EMSchur3D< Eigen::PardisoLDLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Symmetric&Eigen::Lower > >::NewSP();
39
+    //auto EM3D = EMSchur3D< Eigen::PardisoLDLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Symmetric&Eigen::Lower > >::NewSP();
40
 
40
 
41
     // Eigen built-in
41
     // Eigen built-in
42
-    //auto EM3D = EMSchur3D< Eigen::SparseLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::NewSP();
42
+    auto EM3D = EMSchur3D< Eigen::SparseLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::NewSP();
43
 
43
 
44
     // UmfPack
44
     // UmfPack
45
     //auto EM3D = EMSchur3D< Eigen::UmfPackLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::NewSP();
45
     //auto EM3D = EMSchur3D< Eigen::UmfPackLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::NewSP();

+ 9
- 6
include/EMSchur3D.h View File

36
 #include <Eigen/PaStiXSupport>
36
 #include <Eigen/PaStiXSupport>
37
 #endif
37
 #endif
38
 
38
 
39
+#ifdef HAVE_PARDISO
39
 #include <Eigen/PardisoSupport>
40
 #include <Eigen/PardisoSupport>
41
+#endif
40
 
42
 
41
 //#include "CSymSimplicialCholesky.h"
43
 //#include "CSymSimplicialCholesky.h"
42
 
44
 
449
             std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
451
             std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
450
         }
452
         }
451
     }
453
     }
452
-    
453
-// Pardiso 
454
+
455
+#ifdef HAVE_PARDISO
454
     template<>
456
     template<>
455
     void EMSchur3D< Eigen::PardisoLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::BuildCDirectSolver() {
457
     void EMSchur3D< Eigen::PardisoLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::BuildCDirectSolver() {
456
         CSolver = new Eigen::PardisoLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > [Omegas.size()];
458
         CSolver = new Eigen::PardisoLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > [Omegas.size()];
457
         for (int iw=0; iw<Omegas.size(); ++iw) {
459
         for (int iw=0; iw<Omegas.size(); ++iw) {
458
             jsw_timer timer;
460
             jsw_timer timer;
459
             timer.begin();
461
             timer.begin();
460
-        
462
+
461
             //CSolver[iw].pardisoParameterArray()[2] = 28; // OMP_NUM_THREADS?
463
             //CSolver[iw].pardisoParameterArray()[2] = 28; // OMP_NUM_THREADS?
462
-            
464
+
463
             /*  Complex system */
465
             /*  Complex system */
464
             std::cout << "PardisoLU pattern analyzing C_" << iw << ",";
466
             std::cout << "PardisoLU pattern analyzing C_" << iw << ",";
465
             std::cout.flush();
467
             std::cout.flush();
481
         for (int iw=0; iw<Omegas.size(); ++iw) {
483
         for (int iw=0; iw<Omegas.size(); ++iw) {
482
             jsw_timer timer;
484
             jsw_timer timer;
483
             timer.begin();
485
             timer.begin();
484
-        
486
+
485
             //CSolver[iw].pardisoParameterArray()[2] = 28; // OMP_NUM_THREADS?
487
             //CSolver[iw].pardisoParameterArray()[2] = 28; // OMP_NUM_THREADS?
486
             //Eigen::SparseMatrix<Complex> Csym = Cvec[iw].selfadjointView< Eigen::Lower >();
488
             //Eigen::SparseMatrix<Complex> Csym = Cvec[iw].selfadjointView< Eigen::Lower >();
487
-            
489
+
488
             /*  Complex system */
490
             /*  Complex system */
489
             std::cout << "PardisoLDLT pattern analyzing C_" << iw << ",";
491
             std::cout << "PardisoLDLT pattern analyzing C_" << iw << ",";
490
             std::cout.flush();
492
             std::cout.flush();
501
             std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
503
             std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
502
         }
504
         }
503
     }
505
     }
506
+#endif
504
 
507
 
505
 #ifdef HAVE_UMFPACK
508
 #ifdef HAVE_UMFPACK
506
     // Umfpack only seems to work when LOWER and UPPER are set to 1. Workarounds this have not been found.
509
     // Umfpack only seems to work when LOWER and UPPER are set to 1. Workarounds this have not been found.

Loading…
Cancel
Save