Browse Source

Working with Lemma 0.2

master
Trevor Irons 6 years ago
parent
commit
d7b4d91ca5
4 changed files with 69 additions and 54 deletions
  1. 6
    1
      examples/EMSchur3D.cpp
  2. 54
    44
      include/EMSchur3D.h
  3. 4
    4
      include/EMSchur3DBase.h
  4. 5
    5
      src/EMSchur3DBase.cpp

+ 6
- 1
examples/EMSchur3D.cpp View File

29
 int main( int argc, char** argv ) {
29
 int main( int argc, char** argv ) {
30
 
30
 
31
     // BiCGSTAB Diagonal preconditioner
31
     // BiCGSTAB Diagonal preconditioner
32
-    auto EM3D = EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::NewSP();
32
+    auto EM3D = EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::IncompleteLUT<Complex> > >::NewSP();
33
 
33
 
34
     if (argc < 3) {
34
     if (argc < 3) {
35
         std::cout << "EMSchur3D  <rgrid>  <1dmod>  <aemsurvey> " << std::endl;
35
         std::cout << "EMSchur3D  <rgrid>  <1dmod>  <aemsurvey> " << std::endl;
36
         exit( EXIT_SUCCESS );
36
         exit( EXIT_SUCCESS );
37
     }
37
     }
38
 
38
 
39
+    #ifdef LEMMAUSEOMP
40
+    std::cout << "Eigen is using OpenMP!!!\n";
41
+    Eigen::initParallel();
42
+    #endif
43
+
39
     ///////////////////////////////////////////////////
44
     ///////////////////////////////////////////////////
40
     // Read in Grid
45
     // Read in Grid
41
     auto GridRead = RectilinearGridReader::NewSP();
46
     auto GridRead = RectilinearGridReader::NewSP();

+ 54
- 44
include/EMSchur3D.h View File

54
          */
54
          */
55
         static std::shared_ptr< EMSchur3D > NewSP() {
55
         static std::shared_ptr< EMSchur3D > NewSP() {
56
             return std::make_shared< EMSchur3D<Solver> >( ctor_key() );
56
             return std::make_shared< EMSchur3D<Solver> >( ctor_key() );
57
-            //return std::make_shared< EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > > >( ctor_key() ) ;
57
+            //return std::make_shared< EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::RowMajor> > > >( ctor_key() ) ;
58
         }
58
         }
59
 
59
 
60
         /** Default protected constructor, use New */
60
         /** Default protected constructor, use New */
142
     template < class Solver >
142
     template < class Solver >
143
     void EMSchur3D<Solver>::SolveSource ( std::shared_ptr<DipoleSource> Source, const int& isource ) {
143
     void EMSchur3D<Solver>::SolveSource ( std::shared_ptr<DipoleSource> Source, const int& isource ) {
144
 
144
 
145
+        std::cout << "In vanilla SolveSource" << std::endl;
146
+
145
         //  figure out which omega we are working with
147
         //  figure out which omega we are working with
146
         int iw = -1;
148
         int iw = -1;
147
         for (int iiw=0; iiw<Omegas.size(); ++iiw) {
149
         for (int iiw=0; iiw<Omegas.size(); ++iiw) {
168
         FillPoints(dpoint);
170
         FillPoints(dpoint);
169
         PrimaryField(Source, dpoint);
171
         PrimaryField(Source, dpoint);
170
 
172
 
173
+        std::cout << "Done with primary field" << std::endl;
174
+
171
         // Allocate a ton of memory
175
         // Allocate a ton of memory
172
         VectorXcr Phi    = VectorXcr::Zero(uns);
176
         VectorXcr Phi    = VectorXcr::Zero(uns);
173
         VectorXcr ms(unx+uny+unz);  // mu sigma
177
         VectorXcr ms(unx+uny+unz);  // mu sigma
179
         VectorXcr E0     = VectorXcr::Zero(unx+uny+unz);
183
         VectorXcr E0     = VectorXcr::Zero(unx+uny+unz);
180
 
184
 
181
         // Lets get cracking
185
         // Lets get cracking
186
+        std::cout << "Filling source terms" << std::endl;
182
         FillSourceTerms(ms, Se, E0, dpoint, Omegas[iw]);
187
         FillSourceTerms(ms, Se, E0, dpoint, Omegas[iw]);
188
+        std::cout << "Done source terms" << std::endl;
183
 
189
 
184
         /////////////////////////////////////////////////
190
         /////////////////////////////////////////////////
185
         // LOG File
191
         // LOG File
187
         logfile += to_string(isource) + std::string(".log");
193
         logfile += to_string(isource) + std::string(".log");
188
         ofstream logio(logfile.c_str());
194
         ofstream logio(logfile.c_str());
189
 
195
 
190
-        logio << *Source << std::endl;
196
+        std::cout << "just logging" << std::endl;
197
+//        logio << *Source << std::endl;
191
         logio << *Grid << std::endl;
198
         logio << *Grid << std::endl;
192
         logio << *LayModel << std::endl;
199
         logio << *LayModel << std::endl;
200
+        std::cout << "dun logging" << std::endl;
193
 
201
 
194
         // solve for RHS
202
         // solve for RHS
195
         int max_it(nx*ny*nz), iter_done(0);
203
         int max_it(nx*ny*nz), iter_done(0);
334
 
342
 
335
     #ifdef HAVE_SUPERLUMT
343
     #ifdef HAVE_SUPERLUMT
336
     template<>
344
     template<>
337
-    void EMSchur3D< Eigen::SuperLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > >::BuildCDirectSolver() {
345
+    void EMSchur3D< Eigen::SuperLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor> > >::BuildCDirectSolver() {
338
 
346
 
339
-        CSolver = new Eigen::SuperLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > [Omegas.size()];
347
+        CSolver = new Eigen::SuperLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor> > [Omegas.size()];
340
 
348
 
341
         for (int iw=0; iw<Omegas.size(); ++iw) {
349
         for (int iw=0; iw<Omegas.size(); ++iw) {
342
             jsw_timer timer;
350
             jsw_timer timer;
374
     #endif
382
     #endif
375
 
383
 
376
     template<>
384
     template<>
377
-    void EMSchur3D< Eigen::SparseLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::COLAMDOrdering<int> > >::BuildCDirectSolver() {
378
-        CSolver = new Eigen::SparseLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::COLAMDOrdering<int> > [Omegas.size()];
385
+    void EMSchur3D< Eigen::SparseLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::COLAMDOrdering<int> > >::BuildCDirectSolver() {
386
+        CSolver = new Eigen::SparseLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::COLAMDOrdering<int> > [Omegas.size()];
379
         for (int iw=0; iw<Omegas.size(); ++iw) {
387
         for (int iw=0; iw<Omegas.size(); ++iw) {
380
             jsw_timer timer;
388
             jsw_timer timer;
381
             timer.begin();
389
             timer.begin();
399
     }
407
     }
400
 
408
 
401
 //     template<>
409
 //     template<>
402
-//     void EMSchur3D< Eigen::CholmodSupernodalLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
403
-//         CSolver = new Eigen::CholmodSupernodalLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > [Omegas.size()];
410
+//     void EMSchur3D< Eigen::CholmodSupernodalLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
411
+//         CSolver = new Eigen::CholmodSupernodalLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > [Omegas.size()];
404
 //         for (int iw=0; iw<Omegas.size(); ++iw) {
412
 //         for (int iw=0; iw<Omegas.size(); ++iw) {
405
 //             Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
413
 //             Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
406
 //             jsw_timer timer;
414
 //             jsw_timer timer;
420
 //     }
428
 //     }
421
 
429
 
422
 //     template<>
430
 //     template<>
423
-//     void EMSchur3D< Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::NaturalOrdering<int> > > ::BuildCDirectSolver() {
424
-//         CSolver = new Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::NaturalOrdering<int> > [Omegas.size()];
431
+//     void EMSchur3D< Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::NaturalOrdering<int> > > ::BuildCDirectSolver() {
432
+//         CSolver = new Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::NaturalOrdering<int> > [Omegas.size()];
425
 //         for (int iw=0; iw<Omegas.size(); ++iw) {
433
 //         for (int iw=0; iw<Omegas.size(); ++iw) {
426
 //             Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
434
 //             Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
427
 //             jsw_timer timer;
435
 //             jsw_timer timer;
441
 //     }
449
 //     }
442
 //
450
 //
443
 //     template<>
451
 //     template<>
444
-//     void EMSchur3D< Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > > ::BuildCDirectSolver() {
445
-//         CSolver = new Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > [Omegas.size()];
452
+//     void EMSchur3D< Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > > ::BuildCDirectSolver() {
453
+//         CSolver = new Eigen::CSymSimplicialLLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > [Omegas.size()];
446
 //         for (int iw=0; iw<Omegas.size(); ++iw) {
454
 //         for (int iw=0; iw<Omegas.size(); ++iw) {
447
 //             //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
455
 //             //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
448
 //             jsw_timer timer;
456
 //             jsw_timer timer;
462
 //     }
470
 //     }
463
 //
471
 //
464
 //     template<>
472
 //     template<>
465
-//     void EMSchur3D< Eigen::CSymSimplicialLDLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > > ::BuildCDirectSolver() {
466
-//         CSolver = new Eigen::CSymSimplicialLDLT< Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > [Omegas.size()];
473
+//     void EMSchur3D< Eigen::CSymSimplicialLDLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > > ::BuildCDirectSolver() {
474
+//         CSolver = new Eigen::CSymSimplicialLDLT< Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower, Eigen::AMDOrdering<int> > [Omegas.size()];
467
 //         for (int iw=0; iw<Omegas.size(); ++iw) {
475
 //         for (int iw=0; iw<Omegas.size(); ++iw) {
468
 //             Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
476
 //             Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
469
 //             jsw_timer timer;
477
 //             jsw_timer timer;
484
 
492
 
485
 
493
 
486
     template<>
494
     template<>
487
-    void EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::IncompleteLUT<Complex> > > ::BuildCDirectSolver() {
488
-        CSolver = new Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::IncompleteLUT<Complex> > [Omegas.size()];
495
+    void EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::IncompleteLUT<Complex> > > ::BuildCDirectSolver() {
496
+        CSolver = new Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::IncompleteLUT<Complex> > [Omegas.size()];
489
         for (int iw=0; iw<Omegas.size(); ++iw) {
497
         for (int iw=0; iw<Omegas.size(); ++iw) {
490
             Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
498
             Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
499
+            CSolver[iw].preconditioner().setDroptol(1e-3);
500
+            //CSolver[iw].preconditioner().setFillfactor(1e3);
491
             jsw_timer timer;
501
             jsw_timer timer;
492
             timer.begin();
502
             timer.begin();
493
             /*  Complex system */
503
             /*  Complex system */
505
     }
515
     }
506
 
516
 
507
     template<>
517
     template<>
508
-    void EMSchur3D< Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > > ::BuildCDirectSolver() {
509
-        CSolver = new Eigen::BiCGSTAB<Eigen::SparseMatrix<Complex, Eigen::ColMajor> > [Omegas.size()];
518
+    void EMSchur3D< Eigen::BiCGSTAB< Eigen::SparseMatrix<Complex, Eigen::RowMajor> > > ::BuildCDirectSolver() {
519
+        CSolver = new Eigen::BiCGSTAB< Eigen::SparseMatrix<Complex, Eigen::RowMajor> > [Omegas.size()];
510
         for (int iw=0; iw<Omegas.size(); ++iw) {
520
         for (int iw=0; iw<Omegas.size(); ++iw) {
511
             Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
521
             Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
512
             jsw_timer timer;
522
             jsw_timer timer;
525
         }
535
         }
526
     }
536
     }
527
 
537
 
528
-    template<>
529
-    void EMSchur3D<   Eigen::ConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
530
-        CSolver = new Eigen::ConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > [Omegas.size()];
531
-        for (int iw=0; iw<Omegas.size(); ++iw) {
532
-            //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
533
-            jsw_timer timer;
534
-            timer.begin();
535
-            /*  Complex system */
536
-            std::cout << "ConjugateGradient pattern analyzing C_" << iw << ",";
537
-            std::cout.flush();
538
-            CSolver[iw].analyzePattern( Cvec[iw] );
539
-            std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
540
-            // factorize
541
-            timer.begin();
542
-            std::cout << "ConjugateGradient factorising C_" << iw << ", ";
543
-            std::cout.flush();
544
-            CSolver[iw].factorize( Cvec[iw] );
545
-            std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
546
-        }
547
-    }
538
+//     template<>
539
+//     void EMSchur3D<   Eigen::ConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
540
+//         CSolver = new Eigen::ConjugateGradient<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > [Omegas.size()];
541
+//         for (int iw=0; iw<Omegas.size(); ++iw) {
542
+//             //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
543
+//             jsw_timer timer;
544
+//             timer.begin();
545
+//             /*  Complex system */
546
+//             std::cout << "ConjugateGradient pattern analyzing C_" << iw << ",";
547
+//             std::cout.flush();
548
+//             CSolver[iw].analyzePattern( Cvec[iw] );
549
+//             std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
550
+//             // factorize
551
+//             timer.begin();
552
+//             std::cout << "ConjugateGradient factorising C_" << iw << ", ";
553
+//             std::cout.flush();
554
+//             CSolver[iw].factorize( Cvec[iw] );
555
+//             std::cout << " done in " << timer.end() / 60. << " [m]" << std::endl;
556
+//         }
557
+//     }
548
 
558
 
549
 //     template<>
559
 //     template<>
550
-//     void EMSchur3D<   Eigen::PastixLLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
551
-//         CSolver = new Eigen::PastixLLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > [Omegas.size()];
560
+//     void EMSchur3D<   Eigen::PastixLLT<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
561
+//         CSolver = new Eigen::PastixLLT<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > [Omegas.size()];
552
 //         //MPI_Init(NULL, NULL);
562
 //         //MPI_Init(NULL, NULL);
553
 //         for (int iw=0; iw<Omegas.size(); ++iw) {
563
 //         for (int iw=0; iw<Omegas.size(); ++iw) {
554
 //             //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
564
 //             //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
569
 //     }
579
 //     }
570
 //
580
 //
571
 //     template<>
581
 //     template<>
572
-//     void EMSchur3D<   Eigen::PastixLDLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
573
-//         CSolver = new Eigen::PastixLDLT<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, Eigen::Lower > [Omegas.size()];
582
+//     void EMSchur3D<   Eigen::PastixLDLT<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > > ::BuildCDirectSolver() {
583
+//         CSolver = new Eigen::PastixLDLT<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, Eigen::Lower > [Omegas.size()];
574
 //         //MPI_Init(NULL, NULL);
584
 //         //MPI_Init(NULL, NULL);
575
 //         for (int iw=0; iw<Omegas.size(); ++iw) {
585
 //         for (int iw=0; iw<Omegas.size(); ++iw) {
576
 //             //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
586
 //             //Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
592
 //     }
602
 //     }
593
 //
603
 //
594
 //     template<>
604
 //     template<>
595
-//     void EMSchur3D<   Eigen::PastixLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, true > > ::BuildCDirectSolver() {
596
-//         CSolver = new Eigen::PastixLU<Eigen::SparseMatrix<Complex, Eigen::ColMajor>, true > [Omegas.size()];
605
+//     void EMSchur3D<   Eigen::PastixLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, true > > ::BuildCDirectSolver() {
606
+//         CSolver = new Eigen::PastixLU<Eigen::SparseMatrix<Complex, Eigen::RowMajor>, true > [Omegas.size()];
597
 //         //MPI_Init(NULL, NULL);
607
 //         //MPI_Init(NULL, NULL);
598
 //         for (int iw=0; iw<Omegas.size(); ++iw) {
608
 //         for (int iw=0; iw<Omegas.size(); ++iw) {
599
 //             Csym = Cvec[iw].selfadjointView<Eigen::Lower>();
609
 //             Csym = Cvec[iw].selfadjointView<Eigen::Lower>();

+ 4
- 4
include/EMSchur3DBase.h View File

357
 
357
 
358
     /** Matrix objects representing discrete operators
358
     /** Matrix objects representing discrete operators
359
      */
359
      */
360
-    Eigen::SparseMatrix<Complex, Eigen::ColMajor>*                  Cvec;
361
-    Eigen::SparseMatrix<Complex, Eigen::ColMajor>                   D;
360
+    Eigen::SparseMatrix<Complex, Eigen::RowMajor>*                  Cvec;
361
+    Eigen::SparseMatrix<Complex, Eigen::RowMajor>                   D;
362
 
362
 
363
     /** Squeezed matrices for reduced phi grid
363
     /** Squeezed matrices for reduced phi grid
364
      */
364
      */
365
-    Eigen::SparseMatrix<Complex, Eigen::ColMajor>*                  Cvec_s;
366
-    Eigen::SparseMatrix<Complex, Eigen::ColMajor>                   D_s;
365
+    Eigen::SparseMatrix<Complex, Eigen::RowMajor>*                  Cvec_s;
366
+    Eigen::SparseMatrix<Complex, Eigen::RowMajor>                   D_s;
367
 
367
 
368
     /** number of cells in x, set when RectilinearGrid is attached */
368
     /** number of cells in x, set when RectilinearGrid is attached */
369
     int nx;
369
     int nx;

+ 5
- 5
src/EMSchur3DBase.cpp View File

35
     // ====================  FRIEND METHODS  =====================
35
     // ====================  FRIEND METHODS  =====================
36
 
36
 
37
     std::ostream &operator<<(std::ostream &stream, const EMSchur3DBase &ob) {
37
     std::ostream &operator<<(std::ostream &stream, const EMSchur3DBase &ob) {
38
-        stream << ob.Serialize()  << "\n---\n"; // End of doc ---
38
+        stream << ob.Serialize()  << "\n";
39
         return stream;
39
         return stream;
40
     }
40
     }
41
 
41
 
866
         // First set up grid stuff
866
         // First set up grid stuff
867
         if (Cvec) delete [] Cvec;
867
         if (Cvec) delete [] Cvec;
868
         //if (CSolver) delete [] CSolver;
868
         //if (CSolver) delete [] CSolver;
869
-        Cvec    = new Eigen::SparseMatrix<Complex> [Omegas.size()];
869
+        Cvec    = new Eigen::SparseMatrix<Complex, Eigen::RowMajor> [Omegas.size()];
870
 
870
 
871
         //#ifdef LEMMAUSEOMP
871
         //#ifdef LEMMAUSEOMP
872
         //#pragma omp parallel for schedule(static, 1)
872
         //#pragma omp parallel for schedule(static, 1)
1490
         Setup();
1490
         Setup();
1491
 
1491
 
1492
         std::cout << "Entering parallel solve" << std::endl;
1492
         std::cout << "Entering parallel solve" << std::endl;
1493
-        #ifdef LEMMAUSEOMP
1494
-        #pragma omp parallel for schedule(static,1)
1495
-        #endif
1493
+        //#ifdef LEMMAUSEOMP
1494
+        //#pragma omp parallel for schedule(static,1)
1495
+        //#endif
1496
         for (int isource=0; isource<Survey->GetNumberOfSources(); ++isource) {
1496
         for (int isource=0; isource<Survey->GetNumberOfSources(); ++isource) {
1497
             SolveSource(Survey->GetSource(isource), isource);
1497
             SolveSource(Survey->GetSource(isource), isource);
1498
         }
1498
         }

Loading…
Cancel
Save