File: LevelWrap.cpp

package info (click to toggle)
trilinos 12.18.1-2
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 825,604 kB
  • sloc: cpp: 3,352,065; ansic: 432,926; fortran: 169,495; xml: 61,903; python: 40,836; sh: 32,697; makefile: 26,612; javascript: 8,535; perl: 7,136; f90: 6,372; csh: 4,183; objc: 2,620; lex: 1,469; lisp: 810; yacc: 497; awk: 364; ml: 281; php: 145
file content (488 lines) | stat: -rw-r--r-- 21,578 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
// @HEADER
//
// ***********************************************************************
//
//        MueLu: A package for multigrid based preconditioning
//                  Copyright 2012 Sandia Corporation
//
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
// the U.S. Government retains certain rights in this software.
//
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions are
// met:
//
// 1. Redistributions of source code must retain the above copyright
// notice, this list of conditions and the following disclaimer.
//
// 2. Redistributions in binary form must reproduce the above copyright
// notice, this list of conditions and the following disclaimer in the
// documentation and/or other materials provided with the distribution.
//
// 3. Neither the name of the Corporation nor the names of the
// contributors may be used to endorse or promote products derived from
// this software without specific prior written permission.
//
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
//
// Questions? Contact
//                    Jonathan Hu       (jhu@sandia.gov)
//                    Andrey Prokopenko (aprokop@sandia.gov)
//                    Ray Tuminaro      (rstumin@sandia.gov)
//
// ***********************************************************************
//
// @HEADER
#include <iostream>
#include <stdexcept>

// Teuchos
#include <Teuchos_StandardCatchMacros.hpp>

// Galeri
#include <Galeri_XpetraParameters.hpp>
#include <Galeri_XpetraProblemFactory.hpp>
#include <Galeri_XpetraUtils.hpp>
#include <Galeri_XpetraMaps.hpp>

#include <Kokkos_DefaultNode.hpp> // For Epetra only runs this points to FakeKokkos in Xpetra

#include "Xpetra_ConfigDefs.hpp"
#include <Xpetra_Map.hpp>
#include <Xpetra_MultiVector.hpp>
#include <Xpetra_MatrixMatrix.hpp>

#include <MueLu.hpp>
#include <MueLu_Level.hpp>
#include <MueLu_MLParameterListInterpreter.hpp>
#include <MueLu_CreateXpetraPreconditioner.hpp>

// Belos
#ifdef HAVE_MUELU_BELOS
#include "BelosConfigDefs.hpp"
#include "BelosLinearProblem.hpp"
#include "BelosSolverFactory.hpp"
#include "BelosXpetraAdapter.hpp"
#include "BelosMueLuAdapter.hpp"
#endif

using Teuchos::RCP;

//----------------------------------------------------------------------------------------------------------
//
// This example demonstrates how to use MueLu in a fashion that looks like ML's LevelWrap
//
// In this example, we suppose that the user provides a fine matrix A and P & R operators.  These are
// given to MueLu which then forms the next finest matrix and makes a hierarchy from there.
//
//----------------------------------------------------------------------------------------------------------

const std::string thickSeparator = "==========================================================================================================================";
const std::string thinSeparator  = "--------------------------------------------------------------------------------------------------------------------------";

const std::string prefSeparator = "=====================================";

namespace MueLuExamples {

#ifdef HAVE_MUELU_BELOS
  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
  void solve_system_hierarchy(Xpetra::UnderlyingLib& lib,
                              Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>>&   A,
                              Teuchos::RCP<Xpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node>>&   X,
                              Teuchos::RCP<Xpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node>>&   B,
                              Teuchos::RCP<MueLu::Hierarchy<Scalar,LocalOrdinal,GlobalOrdinal,Node>>& H,
                              Teuchos::RCP<Teuchos::ParameterList>& SList) {
#include "MueLu_UseShortNames.hpp"
    using Teuchos::rcp;

    typedef Xpetra::MultiVector <Scalar,LocalOrdinal,GlobalOrdinal,Node> MV;
    typedef Belos::OperatorT<MV>                                         OP;

    // Construct a Belos LinearProblem object
    RCP<OP> belosOp   = rcp(new Belos::XpetraOp<Scalar,LocalOrdinal,GlobalOrdinal,Node>(A));
    RCP<OP> belosPrec = rcp(new Belos::MueLuOp <Scalar,LocalOrdinal,GlobalOrdinal,Node>(H));

    RCP<Belos::LinearProblem<Scalar, MV, OP> > belosProblem =
        rcp(new Belos::LinearProblem<Scalar, MV, OP>(belosOp, X, B));
    belosProblem->setRightPrec(belosPrec);
    belosProblem->setProblem(X,B);
    Belos::SolverFactory<Scalar, MV, OP> BelosFactory;
    RCP<Belos::SolverManager<Scalar, MV, OP> > BelosSolver =
        BelosFactory.create(std::string("CG"), SList);
    BelosSolver->setProblem(belosProblem);
    Belos::ReturnType result = BelosSolver->solve();
    TEUCHOS_TEST_FOR_EXCEPTION(result == Belos::Unconverged, std::runtime_error, "Belos failed to converge");
  }

  // --------------------------------------------------------------------------------------
  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
  void solve_system_list(Xpetra::UnderlyingLib& lib,
                         Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>>& A,
                         Teuchos::RCP<Xpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node>>& X,
                         Teuchos::RCP<Xpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node>>& B,
                         Teuchos::ParameterList& MueLuList,
                         Teuchos::RCP<Teuchos::ParameterList>& SList) {
#include "MueLu_UseShortNames.hpp"
    using Teuchos::rcp;

    if(lib == Xpetra::UseEpetra) {MueLuList.set("use kokkos refactor", false);}
    Teuchos::RCP<MueLu::Hierarchy<Scalar,LocalOrdinal,GlobalOrdinal,Node> > H =
        MueLu::CreateXpetraPreconditioner<Scalar,LocalOrdinal,GlobalOrdinal,Node>(A, MueLuList);

    typedef Xpetra::MultiVector <Scalar,LocalOrdinal,GlobalOrdinal,Node> MV;
    typedef Belos::OperatorT<MV>                                         OP;

    // Construct a Belos LinearProblem object
    RCP<OP> belosOp   = rcp(new Belos::XpetraOp<Scalar,LocalOrdinal,GlobalOrdinal,Node>(A));
    RCP<OP> belosPrec = rcp(new Belos::MueLuOp <Scalar,LocalOrdinal,GlobalOrdinal,Node>(H));

    RCP<Belos::LinearProblem<Scalar, MV, OP> > belosProblem =
        rcp(new Belos::LinearProblem<Scalar, MV, OP>(belosOp, X, B));
    belosProblem->setRightPrec(belosPrec);
    belosProblem->setProblem(X,B);

    Belos::SolverFactory<SC, MV, OP> BelosFactory;
    Teuchos::RCP<Belos::SolverManager<SC, MV, OP> > BelosSolver = BelosFactory.create(std::string("CG"), SList);
    BelosSolver->setProblem(belosProblem);
    Belos::ReturnType result = BelosSolver->solve();
    TEUCHOS_TEST_FOR_EXCEPTION(result == Belos::Unconverged, std::runtime_error, "Belos failed to converge");
  }
#endif


  // --------------------------------------------------------------------------------------
  // This routine generate's the user's original A matrix and nullspace
  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
  void generate_user_matrix_and_nullspace(std::string& matrixType,
                                          Xpetra::UnderlyingLib& lib,
                                          Teuchos::ParameterList& galeriList,
                                          Teuchos::RCP<const Teuchos::Comm<int>>& comm,
                                          Teuchos::RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>>&      A,
                                          Teuchos::RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>>& nullspace) {
#include "MueLu_UseShortNames.hpp"

    RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
    Teuchos::FancyOStream& out = *fancy;

    RCP<const Map>   map;
    RCP<MultiVector> coordinates;
    if (matrixType == "Laplace1D") {
      map = Galeri::Xpetra::CreateMap<LO, GO, Node>(lib, "Cartesian1D", comm, galeriList);
      coordinates = Galeri::Xpetra::Utils::CreateCartesianCoordinates<SC,LO,GO,Map,MultiVector>("1D", map, galeriList);

    } else if (matrixType == "Laplace2D" || matrixType == "Star2D" || matrixType == "BigStar2D" || matrixType == "Elasticity2D") {
      map = Galeri::Xpetra::CreateMap<LO, GO, Node>(lib, "Cartesian2D", comm, galeriList);
      coordinates = Galeri::Xpetra::Utils::CreateCartesianCoordinates<SC,LO,GO,Map,MultiVector>("2D", map, galeriList);

    } else if (matrixType == "Laplace3D" || matrixType == "Brick3D" || matrixType == "Elasticity3D") {
      map = Galeri::Xpetra::CreateMap<LO, GO, Node>(lib, "Cartesian3D", comm, galeriList);
      coordinates = Galeri::Xpetra::Utils::CreateCartesianCoordinates<SC,LO,GO,Map,MultiVector>("3D", map, galeriList);
    }

    // Expand map to do multiple DOF per node for block problems
    if (matrixType == "Elasticity2D" || matrixType == "Elasticity3D")
      map = Xpetra::MapFactory<LO,GO,Node>::Build(map, (matrixType == "Elasticity2D" ? 2 : 3));

    out << "Processor subdomains in x direction: " << galeriList.get<GO>("mx") << std::endl
        << "Processor subdomains in y direction: " << galeriList.get<GO>("my") << std::endl
        << "Processor subdomains in z direction: " << galeriList.get<GO>("mz") << std::endl
        << "========================================================" << std::endl;

    RCP<Galeri::Xpetra::Problem<Map,CrsMatrixWrap,MultiVector> > Pr =
      Galeri::Xpetra::BuildProblem<SC,LO,GO,Map,CrsMatrixWrap,MultiVector>(matrixType, map, galeriList);

    A = Pr->BuildMatrix();

    if (matrixType == "Elasticity2D" || matrixType == "Elasticity3D") {
      nullspace = Pr->BuildNullspace();
      A->SetFixedBlockSize((matrixType == "Elasticity2D") ? 2 : 3);
    }
  }
}


// --------------------------------------------------------------------------------------
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
int main_(Teuchos::CommandLineProcessor &clp, Xpetra::UnderlyingLib lib, int argc, char *argv[]) {
#include <MueLu_UseShortNames.hpp>
  using Teuchos::TimeMonitor;

  bool success = false;
  bool verbose = true;
  try {
#if defined(HAVE_MUELU_SERIAL) && defined(HAVE_TPETRA_INST_INT_INT)
    RCP<const Teuchos::Comm<int> > comm = Teuchos::DefaultComm<int>::getComm();

    RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
    Teuchos::FancyOStream& out = *fancy;

    typedef Teuchos::ScalarTraits<SC> STS;

    // =========================================================================
    // Parameters initialization
    // =========================================================================
    GO nx = 100, ny = 100, nz = 100;
    Galeri::Xpetra::Parameters<GO> galeriParameters(clp, nx, ny, nz, "Laplace2D"); // manage parameters of the test case
    Xpetra::Parameters             xpetraParameters(clp);                          // manage parameters of Xpetra

    switch (clp.parse(argc, argv)) {
      case Teuchos::CommandLineProcessor::PARSE_HELP_PRINTED:        return EXIT_SUCCESS; break;
      case Teuchos::CommandLineProcessor::PARSE_ERROR:
      case Teuchos::CommandLineProcessor::PARSE_UNRECOGNIZED_OPTION: return EXIT_FAILURE; break;
      case Teuchos::CommandLineProcessor::PARSE_SUCCESSFUL:                               break;
    }

    Teuchos::ParameterList galeriList = galeriParameters.GetParameterList();
    out << thickSeparator << std::endl << xpetraParameters << galeriParameters;

    // =========================================================================
    // Problem construction
    // =========================================================================
    RCP<const Map>   map;
    RCP<Matrix> A,P,R, Ac;
    RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal,Node> > nullspace;
    std::string matrixType = galeriParameters.GetMatrixType();
    MueLuExamples::generate_user_matrix_and_nullspace(matrixType,lib,galeriList,comm,A,nullspace);
    map=A->getRowMap();

    // =========================================================================
    // Setups and solves
    // =========================================================================
    RCP<Vector> X = VectorFactory::Build(map);
    RCP<Vector> B = VectorFactory::Build(map);
    B->setSeed(846930886);
    B->randomize();
    RCP<TimeMonitor> tm;

#ifdef HAVE_MUELU_BELOS
    // Belos Options
    RCP<Teuchos::ParameterList> SList = rcp(new Teuchos::ParameterList );
    SList->set("Verbosity",Belos::Errors + Belos::Warnings + Belos::StatusTestDetails);
    SList->set("Output Frequency",10);
    SList->set("Output Style",Belos::Brief);
    SList->set("Maximum Iterations",200);
    SList->set("Convergence Tolerance",1e-10);
#endif


    // =========================================================================
    // Solve #1 (standard MueLu)
    // =========================================================================
    out << thickSeparator << std::endl;
    out << prefSeparator << " Solve 1: Standard "<< prefSeparator <<std::endl;
    {
      // Use an ML-style parameter list for variety
      Teuchos::ParameterList MLList;
      MLList.set("ML output", 10);
      MLList.set("use kokkos refactor", false);
      MLList.set("coarse: type","Amesos-Superlu");
#ifdef HAVE_AMESOS2_KLU2
      MLList.set("coarse: type","Amesos-KLU");
#endif
      MLParameterListInterpreter mueLuFactory(MLList);
      RCP<Hierarchy> H = mueLuFactory.CreateHierarchy();
      Teuchos::RCP<FactoryManagerBase> LevelFactory = mueLuFactory.GetFactoryManager(1);
      H->setlib(lib);
      H->AddNewLevel();
      H->GetLevel(1)->Keep("Nullspace",LevelFactory->GetFactory("Nullspace").get());
      H->GetLevel(0)->Set("A", A);
      mueLuFactory.SetupHierarchy(*H);

#ifdef HAVE_MUELU_BELOS
      // Solve
      MueLuExamples::solve_system_hierarchy<Scalar,LocalOrdinal,GlobalOrdinal,Node>(lib,A,X,B,H,SList);
#endif

      // Extract R,P & Ac for LevelWrap Usage
      H->GetLevel(1)->Get("R",R);
      H->GetLevel(1)->Get("P",P);
      H->GetLevel(1)->Get("A",Ac);

      // extract coarse level null space from level 1 that we have to inject for the next runs...
      nullspace = H->GetLevel(1)->template Get<RCP<MultiVector> >("Nullspace",LevelFactory->GetFactory("Nullspace").get());
    }
    out << thickSeparator << std::endl;

    // =========================================================================
    // Solve #2 (level wrap, the long way, using pre-done Ac)
    // =========================================================================
    out << thickSeparator << std::endl;
    out << prefSeparator << " Solve 2: LevelWrap, Long Way, P, R, Ac "<< prefSeparator <<std::endl;
    {
      // Start w/ an ML-style parameter list
      Teuchos::ParameterList MLList;
      MLList.set("ML output", 10);
      MLList.set("use kokkos refactor", false);
      MLList.set("coarse: type","Amesos-Superlu");
#ifdef HAVE_AMESOS2_KLU2
      MLList.set("coarse: type","Amesos-KLU");
#endif
      FactoryManager M1;
      M1.SetKokkosRefactor(false);
      M1.SetFactory("A",        MueLu::NoFactory::getRCP());
      M1.SetFactory("P",        MueLu::NoFactory::getRCP());
      M1.SetFactory("R",        MueLu::NoFactory::getRCP());

      MLParameterListInterpreter mueLuFactory(MLList);
      mueLuFactory.AddFactoryManager(1, 1, Teuchos::rcpFromRef(M1));
      RCP<Hierarchy> H = mueLuFactory.CreateHierarchy();
      H->setlib(lib);
      H->GetLevel(0)->Set("A", A);
      H->AddNewLevel();
      H->GetLevel(1)->Set("R", R);
      H->GetLevel(1)->Set("P", P);
      H->GetLevel(1)->Set("A", Ac);
      H->GetLevel(1)->Set("Nullspace", nullspace);
      mueLuFactory.SetupHierarchy(*H);

#ifdef HAVE_MUELU_BELOS
      MueLuExamples::solve_system_hierarchy(lib,A,X,B,H,SList);
#endif

    }
    out << thickSeparator << std::endl;


    // =========================================================================
    // Solve #3 (level wrap, the long way, using P, R and nullspace)
    // =========================================================================
    out << thickSeparator << std::endl;
    out << prefSeparator << " Solve 3: LevelWrap, Long Way, P, R "<< prefSeparator <<std::endl;
    {

      // Start w/ an ML-style parameter list
      Teuchos::ParameterList MLList;
      MLList.set("ML output", 10);
      MLList.set("use kokkos refactor", false);
      MLList.set("coarse: type","Amesos-Superlu");
#ifdef HAVE_AMESOS2_KLU2
      MLList.set("coarse: type","Amesos-KLU");
#endif
      FactoryManager M1;
      M1.SetKokkosRefactor(false);
      M1.SetFactory("P",        MueLu::NoFactory::getRCP());
      M1.SetFactory("R",        MueLu::NoFactory::getRCP());

      MLParameterListInterpreter mueLuFactory(MLList);
      mueLuFactory.AddFactoryManager(1, 1, Teuchos::rcpFromRef(M1));
      RCP<Hierarchy> H = mueLuFactory.CreateHierarchy();
      H->setlib(lib);
      H->GetLevel(0)->Set("A", A);
      H->AddNewLevel();
      H->GetLevel(1)->Set("R", R);
      H->GetLevel(1)->Set("P", P);
      H->GetLevel(1)->Set("Nullspace", nullspace);
      mueLuFactory.SetupHierarchy(*H);
#ifdef HAVE_MUELU_BELOS
      MueLuExamples::solve_system_hierarchy(lib,A,X,B,H,SList);
#endif

    }
    out << thickSeparator << std::endl;

    // =========================================================================
    // Solve #4 (level wrap, the fast way, everything)
    // =========================================================================
    out << thickSeparator << std::endl;
    out << prefSeparator << " Solve 4: LevelWrap, Fast Way, P, R, Ac "<< prefSeparator <<std::endl;
    {
      Teuchos::ParameterList MueLuList, level1;
      level1.set("A",Ac);
      level1.set("R",R);
      level1.set("P",P);
      level1.set("Nullspace",nullspace);
      MueLuList.set("level 1",level1);
      MueLuList.set("verbosity","high");
      MueLuList.set("coarse: max size",100);
      MueLuList.set("max levels",4);
#ifdef HAVE_MUELU_BELOS
      MueLuExamples::solve_system_list(lib,A,X,B,MueLuList,SList);
#endif
    }

    // =========================================================================
    // Solve #5 (level wrap, the fast way, P, R + Nullspace)
    // =========================================================================
    out << thickSeparator << std::endl;
    out << prefSeparator << " Solve 5: LevelWrap, Fast Way, P, R "<< prefSeparator <<std::endl;
    {
      Teuchos::ParameterList MueLuList, level1;
      level1.set("R",R);
      level1.set("P",P);
      level1.set("Nullspace",nullspace);
      MueLuList.set("level 1",level1);
      MueLuList.set("verbosity","high");
      MueLuList.set("coarse: max size",100);
#ifdef HAVE_MUELU_BELOS
      MueLuExamples::solve_system_list(lib,A,X,B,MueLuList,SList);
#endif
    }


    // =========================================================================
    // Solve #6 (level wrap, the fast way, P only, explicit transpose)
    // =========================================================================
    out << thickSeparator << std::endl;
    out << prefSeparator << " Solve 6: LevelWrap, Fast Way, P only, explicit transpose "<< prefSeparator <<std::endl;
    {
      Teuchos::ParameterList MueLuList, level1;
      level1.set("P",P);
      level1.set("Nullspace",nullspace);
      MueLuList.set("level 1",level1);
      MueLuList.set("verbosity","high");
      MueLuList.set("coarse: max size",100);
      MueLuList.set("transpose: use implicit",false);
      MueLuList.set("max levels",4);
#ifdef HAVE_MUELU_BELOS
      MueLuExamples::solve_system_list(lib,A,X,B,MueLuList,SList);
#endif
    }


    // =========================================================================
    // Solve #7 (level wrap, the fast way, P only, implicit transpose)
    // =========================================================================
    out << thickSeparator << std::endl;
    out << prefSeparator << " Solve 7: LevelWrap, Fast Way, P only, implicit transpose "<< prefSeparator <<std::endl;
    {
      Teuchos::ParameterList MueLuList, level1;
      level1.set("P",P);
      level1.set("Nullspace",nullspace);
      MueLuList.set("level 1",level1);
      MueLuList.set("verbosity","high");
      MueLuList.set("coarse: max size",100);
      MueLuList.set("transpose: use implicit",true);
      MueLuList.set("max levels",2);
#ifdef HAVE_MUELU_BELOS
      MueLuExamples::solve_system_list(lib,A,X,B,MueLuList,SList);
#endif
    }
#endif
    success = true;
  }
  TEUCHOS_STANDARD_CATCH_STATEMENTS(verbose, std::cerr, success);

  return ( success ? EXIT_SUCCESS : EXIT_FAILURE );
}



//- -- --------------------------------------------------------
#define MUELU_AUTOMATIC_TEST_ETI_NAME main_
#include "MueLu_Test_ETI.hpp"

int main(int argc, char *argv[]) {
  return Automatic_Test_ETI(argc,argv);
}