MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_Ifpack2Smoother_def.hpp
Go to the documentation of this file.
1// @HEADER
2//
3// ***********************************************************************
4//
5// MueLu: A package for multigrid based preconditioning
6// Copyright 2012 Sandia Corporation
7//
8// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9// the U.S. Government retains certain rights in this software.
10//
11// Redistribution and use in source and binary forms, with or without
12// modification, are permitted provided that the following conditions are
13// met:
14//
15// 1. Redistributions of source code must retain the above copyright
16// notice, this list of conditions and the following disclaimer.
17//
18// 2. Redistributions in binary form must reproduce the above copyright
19// notice, this list of conditions and the following disclaimer in the
20// documentation and/or other materials provided with the distribution.
21//
22// 3. Neither the name of the Corporation nor the names of the
23// contributors may be used to endorse or promote products derived from
24// this software without specific prior written permission.
25//
26// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37//
38// Questions? Contact
39// Jonathan Hu (jhu@sandia.gov)
40// Andrey Prokopenko (aprokop@sandia.gov)
41// Ray Tuminaro (rstumin@sandia.gov)
42//
43// ***********************************************************************
44//
45// @HEADER
46#ifndef MUELU_IFPACK2SMOOTHER_DEF_HPP
47#define MUELU_IFPACK2SMOOTHER_DEF_HPP
48
49#include "MueLu_ConfigDefs.hpp"
50
51#if defined(HAVE_MUELU_IFPACK2)
52
53#include <Teuchos_ParameterList.hpp>
54
55#include <Tpetra_RowMatrix.hpp>
56
57#include <Ifpack2_Chebyshev.hpp>
58#include <Ifpack2_RILUK.hpp>
59#include <Ifpack2_Relaxation.hpp>
60#include <Ifpack2_ILUT.hpp>
61#include <Ifpack2_BlockRelaxation.hpp>
62#include <Ifpack2_Factory.hpp>
63#include <Ifpack2_Parameters.hpp>
64
65#include <Xpetra_BlockedCrsMatrix.hpp>
66#include <Xpetra_CrsMatrix.hpp>
67#include <Xpetra_CrsMatrixWrap.hpp>
68#include <Xpetra_TpetraBlockCrsMatrix.hpp>
69#include <Xpetra_Matrix.hpp>
70#include <Xpetra_MatrixMatrix.hpp>
71#include <Xpetra_MultiVectorFactory.hpp>
72#include <Xpetra_TpetraMultiVector.hpp>
73
74#include <Tpetra_BlockCrsMatrix_Helpers.hpp>
75
77#include "MueLu_Level.hpp"
78#include "MueLu_Utilities.hpp"
79#include "MueLu_Monitor.hpp"
80#include "MueLu_Aggregates.hpp"
81
82
83#ifdef HAVE_MUELU_INTREPID2
86#include "Intrepid2_Basis.hpp"
87#include "Kokkos_DynRankView.hpp"
88#endif
89
90
91namespace MueLu {
92
93 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
94 Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Ifpack2Smoother(const std::string& type, const Teuchos::ParameterList& paramList, const LO& overlap)
95 : type_(type), overlap_(overlap)
96 {
97 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
98 bool isSupported = Ifpack2::Factory::isSupported<tRowMatrix>(type_) || (type_ == "LINESMOOTHING_TRIDI_RELAXATION" ||
99 type_ == "LINESMOOTHING_TRIDI RELAXATION" ||
100 type_ == "LINESMOOTHING_TRIDIRELAXATION" ||
101 type_ == "LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
102 type_ == "LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
103 type_ == "LINESMOOTHING_TRIDIAGONALRELAXATION" ||
104 type_ == "LINESMOOTHING_BANDED_RELAXATION" ||
105 type_ == "LINESMOOTHING_BANDED RELAXATION" ||
106 type_ == "LINESMOOTHING_BANDEDRELAXATION" ||
107 type_ == "LINESMOOTHING_BLOCK_RELAXATION" ||
108 type_ == "LINESMOOTHING_BLOCK RELAXATION" ||
109 type_ == "LINESMOOTHING_BLOCKRELAXATION" ||
110 type_ == "TOPOLOGICAL" ||
111 type_ == "AGGREGATE");
112 this->declareConstructionOutcome(!isSupported, "Ifpack2 does not provide the smoother '" + type_ + "'.");
113 if (isSupported)
114 SetParameterList(paramList);
115 }
116
117 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
119 Factory::SetParameterList(paramList);
120
122 // It might be invalid to change parameters after the setup, but it depends entirely on Ifpack implementation.
123 // TODO: I don't know if Ifpack returns an error code or exception or ignore parameters modification in this case...
125 }
126 }
127
128 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
130 std::string prefix = this->ShortClassName() + ": SetPrecParameters";
131 RCP<TimeMonitor> tM = rcp(new TimeMonitor(*this, prefix, Timings0));
132 ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
133
134 paramList.setParameters(list);
135
136 RCP<ParameterList> precList = this->RemoveFactoriesFromList(this->GetParameterList());
137
138 // Do we want an Ifpack2 apply timer?
139 precList->set("timer for apply", this->IsPrint(Timings));
140
141 if(!precList.is_null() && precList->isParameter("partitioner: type") && precList->get<std::string>("partitioner: type") == "linear" &&
142 !precList->isParameter("partitioner: local parts")) {
143 LO matrixBlockSize = 1;
144 int lclSize = A_->getRangeMap()->getLocalNumElements();
145 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
146 if (!matA.is_null()) {
147 lclSize = matA->getLocalNumRows();
148 matrixBlockSize = matA->GetFixedBlockSize();
149 }
150 precList->set("partitioner: local parts", lclSize / matrixBlockSize);
151 }
152
153 prec_->setParameters(*precList);
154
155 paramList.setParameters(*precList); // what about that??
156 }
157
158 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
160 this->Input(currentLevel, "A");
161
162 if (type_ == "LINESMOOTHING_TRIDI_RELAXATION" ||
163 type_ == "LINESMOOTHING_TRIDI RELAXATION" ||
164 type_ == "LINESMOOTHING_TRIDIRELAXATION" ||
165 type_ == "LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
166 type_ == "LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
167 type_ == "LINESMOOTHING_TRIDIAGONALRELAXATION" ||
168 type_ == "LINESMOOTHING_BANDED_RELAXATION" ||
169 type_ == "LINESMOOTHING_BANDED RELAXATION" ||
170 type_ == "LINESMOOTHING_BANDEDRELAXATION" ||
171 type_ == "LINESMOOTHING_BLOCK_RELAXATION" ||
172 type_ == "LINESMOOTHING_BLOCK RELAXATION" ||
173 type_ == "LINESMOOTHING_BLOCKRELAXATION") {
174 this->Input(currentLevel, "CoarseNumZLayers"); // necessary for fallback criterion
175 this->Input(currentLevel, "LineDetection_VertLineIds"); // necessary to feed block smoother
176 }
177 else if (type_ == "BLOCK RELAXATION" ||
178 type_ == "BLOCK_RELAXATION" ||
179 type_ == "BLOCKRELAXATION" ||
180 // Banded
181 type_ == "BANDED_RELAXATION" ||
182 type_ == "BANDED RELAXATION" ||
183 type_ == "BANDEDRELAXATION" ||
184 // Tridiagonal
185 type_ == "TRIDI_RELAXATION" ||
186 type_ == "TRIDI RELAXATION" ||
187 type_ == "TRIDIRELAXATION" ||
188 type_ == "TRIDIAGONAL_RELAXATION" ||
189 type_ == "TRIDIAGONAL RELAXATION" ||
190 type_ == "TRIDIAGONALRELAXATION")
191 {
192 //We need to check for the "partitioner type" = "line"
193 ParameterList precList = this->GetParameterList();
194 if(precList.isParameter("partitioner: type") &&
195 precList.get<std::string>("partitioner: type") == "line") {
196 this->Input(currentLevel, "Coordinates");
197 }
198 }
199 else if (type_ == "TOPOLOGICAL")
200 {
201 // for the topological smoother, we require an element to node map:
202 this->Input(currentLevel, "pcoarsen: element to node map");
203 }
204 else if (type_ == "AGGREGATE")
205 {
206 // Aggregate smoothing needs aggregates
207 this->Input(currentLevel,"Aggregates");
208 }
209 else if (type_ == "HIPTMAIR") {
210 // Hiptmair needs D0 and NodeMatrix
211 this->Input(currentLevel,"NodeMatrix");
212 this->Input(currentLevel,"D0");
214
217 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
219 FactoryMonitor m(*this, "Setup Smoother", currentLevel);
220 A_ = Factory::Get< RCP<Operator> >(currentLevel, "A");
221 ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
222
223 // If the user asked us to convert the matrix into BlockCrsMatrix form, we do that now
224
225 if(paramList.isParameter("smoother: use blockcrsmatrix storage") && paramList.get<bool>("smoother: use blockcrsmatrix storage")) {
226 int blocksize = 1;
227 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
228 if (!matA.is_null())
229 blocksize = matA->GetFixedBlockSize();
230 if (blocksize) {
231 // NOTE: Don't think you can move this out of the if block. You can't. The test MueLu_MeshTyingBlocked_SimpleSmoother_2dof_medium_MPI_1 will fail
232
233 RCP<CrsMatrixWrap> AcrsWrap = rcp_dynamic_cast<CrsMatrixWrap>(A_);
234 if(AcrsWrap.is_null())
235 throw std::runtime_error("Ifpack2Smoother: Cannot convert matrix A to CrsMatrixWrap object.");
236 RCP<CrsMatrix> Acrs = AcrsWrap->getCrsMatrix();
237 if(Acrs.is_null())
238 throw std::runtime_error("Ifpack2Smoother: Cannot extract CrsMatrix from matrix A.");
239 RCP<TpetraCrsMatrix> At = rcp_dynamic_cast<TpetraCrsMatrix>(Acrs);
240 if(At.is_null()) {
241 if(!Xpetra::Helpers<Scalar,LO,GO,Node>::isTpetraBlockCrs(matA))
242 throw std::runtime_error("Ifpack2Smoother: Cannot extract CrsMatrix or BlockCrsMatrix from matrix A.");
243 this->GetOStream(Statistics0) << "Ifpack2Smoother: Using (native) BlockCrsMatrix storage with blocksize "<<blocksize<<std::endl;
244 paramList.remove("smoother: use blockcrsmatrix storage");
245 }
246 else {
247 RCP<Tpetra::BlockCrsMatrix<Scalar, LO, GO, Node> > blockCrs = Tpetra::convertToBlockCrsMatrix(*At->getTpetra_CrsMatrix(),blocksize);
248 RCP<CrsMatrix> blockCrs_as_crs = rcp(new TpetraBlockCrsMatrix(blockCrs));
249 RCP<CrsMatrixWrap> blockWrap = rcp(new CrsMatrixWrap(blockCrs_as_crs));
250 A_ = blockWrap;
251 this->GetOStream(Statistics0) << "Ifpack2Smoother: Using BlockCrsMatrix storage with blocksize "<<blocksize<<std::endl;
252
253 paramList.remove("smoother: use blockcrsmatrix storage");
254 }
255 }
256 }
257
258 if (type_ == "SCHWARZ")
259 SetupSchwarz(currentLevel);
260
261 else if (type_ == "LINESMOOTHING_TRIDI_RELAXATION" ||
262 type_ == "LINESMOOTHING_TRIDI RELAXATION" ||
263 type_ == "LINESMOOTHING_TRIDIRELAXATION" ||
264 type_ == "LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
265 type_ == "LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
266 type_ == "LINESMOOTHING_TRIDIAGONALRELAXATION" ||
267 type_ == "LINESMOOTHING_BANDED_RELAXATION" ||
268 type_ == "LINESMOOTHING_BANDED RELAXATION" ||
269 type_ == "LINESMOOTHING_BANDEDRELAXATION" ||
270 type_ == "LINESMOOTHING_BLOCK_RELAXATION" ||
271 type_ == "LINESMOOTHING_BLOCK RELAXATION" ||
272 type_ == "LINESMOOTHING_BLOCKRELAXATION")
273 SetupLineSmoothing(currentLevel);
274
275 else if (type_ == "BLOCK_RELAXATION" ||
276 type_ == "BLOCK RELAXATION" ||
277 type_ == "BLOCKRELAXATION" ||
278 // Banded
279 type_ == "BANDED_RELAXATION" ||
280 type_ == "BANDED RELAXATION" ||
281 type_ == "BANDEDRELAXATION" ||
282 // Tridiagonal
283 type_ == "TRIDI_RELAXATION" ||
284 type_ == "TRIDI RELAXATION" ||
285 type_ == "TRIDIRELAXATION" ||
286 type_ == "TRIDIAGONAL_RELAXATION" ||
287 type_ == "TRIDIAGONAL RELAXATION" ||
288 type_ == "TRIDIAGONALRELAXATION")
289 SetupBlockRelaxation(currentLevel);
290
291 else if (type_ == "CHEBYSHEV")
292 SetupChebyshev(currentLevel);
293
294 else if (type_ == "TOPOLOGICAL")
295 {
296#ifdef HAVE_MUELU_INTREPID2
297 SetupTopological(currentLevel);
298#else
299 TEUCHOS_TEST_FOR_EXCEPTION(true, std::invalid_argument, "'TOPOLOGICAL' smoother choice requires Intrepid2");
300#endif
301 }
302 else if (type_ == "AGGREGATE")
303 SetupAggregate(currentLevel);
304
305 else if (type_ == "HIPTMAIR")
306 SetupHiptmair(currentLevel);
307
308 else
309 {
310 SetupGeneric(currentLevel);
311 }
312
314
315 this->GetOStream(Statistics1) << description() << std::endl;
316 }
317
318 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
320 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
321
322 bool reusePreconditioner = false;
323 if (this->IsSetup() == true) {
324 // Reuse the constructed preconditioner
325 this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupSchwarz(): Setup() has already been called, assuming reuse" << std::endl;
326
327 bool isTRowMatrix = true;
328 RCP<const tRowMatrix> tA;
329 try {
331 } catch (Exceptions::BadCast&) {
332 isTRowMatrix = false;
333 }
334
335 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(prec_);
336 if (!prec.is_null() && isTRowMatrix) {
337 prec->setMatrix(tA);
338 reusePreconditioner = true;
339 } else {
340 this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupSchwarz(): reuse of this type is not available "
341 "(either failed cast to CanChangeMatrix, or to Tpetra Row Matrix), reverting to full construction" << std::endl;
342 }
343 }
344
345 if (!reusePreconditioner) {
346 ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
347
348 bool isBlockedMatrix = false;
349 RCP<Matrix> merged2Mat;
350
351 std::string sublistName = "subdomain solver parameters";
352 if (paramList.isSublist(sublistName)) {
353 // If we are doing "user" partitioning, we assume that what the user
354 // really wants to do is make tiny little subdomains with one row
355 // assigned to each subdomain. The rows used for these little
356 // subdomains correspond to those in the 2nd block row. Then,
357 // if we overlap these mini-subdomains, we will do something that
358 // looks like Vanka (grabbing all velocities associated with each
359 // each pressure unknown). In addition, we put all Dirichlet points
360 // as a little mini-domain.
361 ParameterList& subList = paramList.sublist(sublistName);
362
363 std::string partName = "partitioner: type";
364 // Pretty sure no one has been using this. Unfortunately, old if
365 // statement (which checked for equality with "user") prevented
366 // anyone from employing other types of Ifpack2 user partition
367 // options. Leaving this and switching if to "vanka user" just
368 // in case some day someone might want to use this.
369 if (subList.isParameter(partName) && subList.get<std::string>(partName) == "vanka user") {
370 isBlockedMatrix = true;
371
372 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
373 TEUCHOS_TEST_FOR_EXCEPTION(bA.is_null(), Exceptions::BadCast,
374 "Matrix A must be of type BlockedCrsMatrix.");
375
376 size_t numVels = bA->getMatrix(0,0)->getLocalNumRows();
377 size_t numPres = bA->getMatrix(1,0)->getLocalNumRows();
378 size_t numRows = rcp_dynamic_cast<Matrix>(A_, true)->getLocalNumRows();
379
380 ArrayRCP<LocalOrdinal> blockSeeds(numRows, Teuchos::OrdinalTraits<LocalOrdinal>::invalid());
381
382 size_t numBlocks = 0;
383 for (size_t rowOfB = numVels; rowOfB < numVels+numPres; ++rowOfB)
384 blockSeeds[rowOfB] = numBlocks++;
385
386 RCP<BlockedCrsMatrix> bA2 = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
387 TEUCHOS_TEST_FOR_EXCEPTION(bA2.is_null(), Exceptions::BadCast,
388 "Matrix A must be of type BlockedCrsMatrix.");
389
390 merged2Mat = bA2->Merge();
391
392 // Add Dirichlet rows to the list of seeds
393 ArrayRCP<const bool> boundaryNodes = Utilities::DetectDirichletRows(*merged2Mat, 0.0);
394 bool haveBoundary = false;
395 for (LO i = 0; i < boundaryNodes.size(); i++)
396 if (boundaryNodes[i]) {
397 // FIXME:
398 // 1. would not this [] overlap with some in the previos blockSeed loop?
399 // 2. do we need to distinguish between pressure and velocity Dirichlet b.c.
400 blockSeeds[i] = numBlocks;
401 haveBoundary = true;
402 }
403 if (haveBoundary)
404 numBlocks++;
405
406 subList.set("partitioner: type","user");
407 subList.set("partitioner: map", blockSeeds);
408 subList.set("partitioner: local parts", as<int>(numBlocks));
409
410 } else {
411 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
412 if (!bA.is_null()) {
413 isBlockedMatrix = true;
414 merged2Mat = bA->Merge();
415 }
416 }
417 }
418
419 RCP<const tRowMatrix> tA;
420 if (isBlockedMatrix == true) tA = Utilities::Op2NonConstTpetraRow(merged2Mat);
422
425
426 prec_->initialize();
427 }
428
429 prec_->compute();
430 }
431
432
433
434 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
436
437 ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
438
439 if (this->IsSetup() == true) {
440 this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupAggregate(): Setup() has already been called" << std::endl;
441 this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupAggregate(): reuse of this type is not available, reverting to full construction" << std::endl;
442 }
443
444 this->GetOStream(Statistics0) << "Ifpack2Smoother: Using Aggregate Smoothing"<<std::endl;
445
446 RCP<Aggregates> aggregates = Factory::Get<RCP<Aggregates> >(currentLevel,"Aggregates");
447
448 RCP<const LOMultiVector> vertex2AggId = aggregates->GetVertex2AggId();
449 ArrayRCP<LO> aggregate_ids = rcp_const_cast<LOMultiVector>(vertex2AggId)->getDataNonConst(0);
450 ArrayRCP<LO> dof_ids;
451
452 // We need to unamalgamate, if the FixedBlockSize > 1
453 LO blocksize = 1;
454 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
455 if (!matA.is_null())
456 blocksize = matA->GetFixedBlockSize();
457 if(blocksize > 1) {
458 dof_ids.resize(aggregate_ids.size() * blocksize);
459 for(LO i=0; i<(LO)aggregate_ids.size(); i++) {
460 for(LO j=0; j<(LO)blocksize; j++)
461 dof_ids[i*blocksize+j] = aggregate_ids[i];
462 }
463 }
464 else {
465 dof_ids = aggregate_ids;
466 }
467
468
469 paramList.set("partitioner: map", dof_ids);
470 paramList.set("partitioner: type", "user");
471 paramList.set("partitioner: overlap", 0);
472 paramList.set("partitioner: local parts", (int)aggregates->GetNumAggregates());
473
474 RCP<const Tpetra::RowMatrix<SC, LO, GO, NO> > tA = Utilities::Op2NonConstTpetraRow(A_);
475
476 type_ = "BLOCKRELAXATION";
479 prec_->initialize();
480 prec_->compute();
481 }
482
483#ifdef HAVE_MUELU_INTREPID2
484 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
486 /*
487
488 basic notion:
489
490 Look for user input indicating topo dimension, something like "topological domain type: {node|edge|face}"
491 Call something like what you can find in Poisson example line 1180 to set seeds for a smoother
492
493 */
494 if (this->IsSetup() == true) {
495 this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupTopological(): Setup() has already been called" << std::endl;
496 this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupTopological(): reuse of this type is not available, reverting to full construction" << std::endl;
497 }
498
499 ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
500
501 typedef typename Node::device_type::execution_space ES;
502
503 typedef Kokkos::DynRankView<LocalOrdinal,typename Node::device_type> FCO; //
504
505 LocalOrdinal lo_invalid = Teuchos::OrdinalTraits<LO>::invalid();
506
507 using namespace std;
508
509 const Teuchos::RCP<FCO> elemToNode = Factory::Get<Teuchos::RCP<FCO> >(currentLevel,"pcoarsen: element to node map");
510
511 string basisString = paramList.get<string>("pcoarsen: hi basis");
512 int degree;
513 // NOTE: To make sure Stokhos works we only instantiate these guys with double. There's a lot
514 // of stuff in the guts of Intrepid2 that doesn't play well with Stokhos as of yet. Here, we only
515 // care about the assignment of basis ordinals to topological entities, so this code is actually
516 // independent of the Scalar type--hard-coding double here won't hurt us.
517 auto basis = MueLuIntrepid::BasisFactory<double,ES>(basisString, degree);
518
519 string topologyTypeString = paramList.get<string>("smoother: neighborhood type");
520 int dimension;
521 if (topologyTypeString == "node")
522 dimension = 0;
523 else if (topologyTypeString == "edge")
524 dimension = 1;
525 else if (topologyTypeString == "face")
526 dimension = 2;
527 else if (topologyTypeString == "cell")
528 dimension = basis->getBaseCellTopology().getDimension();
529 else
530 TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Unrecognized smoother neighborhood type. Supported types are node, edge, face.");
531 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_, true);
532 vector<vector<LocalOrdinal>> seeds;
533 MueLuIntrepid::FindGeometricSeedOrdinals(basis, *elemToNode, seeds, *matA->getRowMap(), *matA->getColMap());
534
535 // Ifpack2 wants the seeds in an array of the same length as the number of local elements,
536 // with local partition #s marked for the ones that are seeds, and invalid for the rest
537 int myNodeCount = matA->getRowMap()->getLocalNumElements();
538 ArrayRCP<LocalOrdinal> nodeSeeds(myNodeCount,lo_invalid);
539 int localPartitionNumber = 0;
540 for (LocalOrdinal seed : seeds[dimension])
541 {
542 nodeSeeds[seed] = localPartitionNumber++;
543 }
544
545 paramList.remove("smoother: neighborhood type");
546 paramList.remove("pcoarsen: hi basis");
547
548 paramList.set("partitioner: map", nodeSeeds);
549 paramList.set("partitioner: type", "user");
550 paramList.set("partitioner: overlap", 1);
551 paramList.set("partitioner: local parts", int(seeds[dimension].size()));
552
553 RCP<const Tpetra::RowMatrix<SC, LO, GO, NO> > tA = Utilities::Op2NonConstTpetraRow(A_);
554
555 type_ = "BLOCKRELAXATION";
556 prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
557 SetPrecParameters();
558 prec_->initialize();
559 prec_->compute();
560 }
561#endif
562
563 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
565 if (this->IsSetup() == true) {
566 this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupLineSmoothing(): Setup() has already been called" << std::endl;
567 this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupLineSmoothing(): reuse of this type is not available, reverting to full construction" << std::endl;
568 }
569
570 ParameterList& myparamList = const_cast<ParameterList&>(this->GetParameterList());
571
572 LO CoarseNumZLayers = Factory::Get<LO>(currentLevel,"CoarseNumZLayers");
573 if (CoarseNumZLayers > 0) {
574 Teuchos::ArrayRCP<LO> TVertLineIdSmoo = Factory::Get< Teuchos::ArrayRCP<LO> >(currentLevel, "LineDetection_VertLineIds");
575
576 // determine number of local parts
577 LO maxPart = 0;
578 for(size_t k = 0; k < Teuchos::as<size_t>(TVertLineIdSmoo.size()); k++) {
579 if(maxPart < TVertLineIdSmoo[k]) maxPart = TVertLineIdSmoo[k];
580 }
581 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_, true);
582 size_t numLocalRows = matA->getLocalNumRows();
583
584 TEUCHOS_TEST_FOR_EXCEPTION(numLocalRows % TVertLineIdSmoo.size() != 0, Exceptions::RuntimeError,
585 "MueLu::Ifpack2Smoother::Setup(): the number of local nodes is incompatible with the TVertLineIdsSmoo.");
586
587 //actualDofsPerNode is the actual number of matrix rows per mesh element.
588 //It is encoded in either the MueLu Level, or in the Xpetra matrix block size.
589 //This value is needed by Ifpack2 to do decoupled block relaxation.
590 int actualDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
591 LO matrixBlockSize = matA->GetFixedBlockSize();
592 if(matrixBlockSize > 1 && actualDofsPerNode > 1)
593 {
594 TEUCHOS_TEST_FOR_EXCEPTION(actualDofsPerNode != matrixBlockSize, Exceptions::RuntimeError,
595 "MueLu::Ifpack2Smoother::Setup(): A is a block matrix but its block size and DOFs/node from partitioner disagree");
596 }
597 else if(matrixBlockSize > 1)
598 {
599 actualDofsPerNode = matrixBlockSize;
600 }
601 myparamList.set("partitioner: PDE equations", actualDofsPerNode);
602
603 if (numLocalRows == Teuchos::as<size_t>(TVertLineIdSmoo.size())) {
604 myparamList.set("partitioner: type","user");
605 myparamList.set("partitioner: map",TVertLineIdSmoo);
606 myparamList.set("partitioner: local parts",maxPart+1);
607 } else {
608 // we assume a constant number of DOFs per node
609 size_t numDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
610
611 // Create a new Teuchos::ArrayRCP<LO> of size numLocalRows and fill it with the corresponding information
612 Teuchos::ArrayRCP<LO> partitionerMap(numLocalRows, Teuchos::OrdinalTraits<LocalOrdinal>::invalid());
613 for (size_t blockRow = 0; blockRow < Teuchos::as<size_t>(TVertLineIdSmoo.size()); ++blockRow)
614 for (size_t dof = 0; dof < numDofsPerNode; dof++)
615 partitionerMap[blockRow * numDofsPerNode + dof] = TVertLineIdSmoo[blockRow];
616 myparamList.set("partitioner: type","user");
617 myparamList.set("partitioner: map",partitionerMap);
618 myparamList.set("partitioner: local parts",maxPart + 1);
619 }
620
621 if (type_ == "LINESMOOTHING_BANDED_RELAXATION" ||
622 type_ == "LINESMOOTHING_BANDED RELAXATION" ||
623 type_ == "LINESMOOTHING_BANDEDRELAXATION")
624 type_ = "BANDEDRELAXATION";
625 else if (type_ == "LINESMOOTHING_TRIDI_RELAXATION" ||
626 type_ == "LINESMOOTHING_TRIDI RELAXATION" ||
627 type_ == "LINESMOOTHING_TRIDIRELAXATION" ||
628 type_ == "LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
629 type_ == "LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
630 type_ == "LINESMOOTHING_TRIDIAGONALRELAXATION")
631 type_ = "TRIDIAGONALRELAXATION";
632 else
633 type_ = "BLOCKRELAXATION";
634 } else {
635 // line detection failed -> fallback to point-wise relaxation
636 this->GetOStream(Runtime0) << "Line detection failed: fall back to point-wise relaxation" << std::endl;
637 myparamList.remove("partitioner: type",false);
638 myparamList.remove("partitioner: map", false);
639 myparamList.remove("partitioner: local parts",false);
640 type_ = "RELAXATION";
641 }
642
643 RCP<const Tpetra::RowMatrix<SC, LO, GO, NO> > tA = Utilities::Op2NonConstTpetraRow(A_);
644
647 prec_->initialize();
648 prec_->compute();
649 }
650
651 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
653 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
654
655 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
656 if (!bA.is_null())
657 A_ = bA->Merge();
658
659 RCP<const tRowMatrix> tA = Utilities::Op2NonConstTpetraRow(A_);
660
661 bool reusePreconditioner = false;
662 if (this->IsSetup() == true) {
663 // Reuse the constructed preconditioner
664 this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupBlockRelaxation(): Setup() has already been called, assuming reuse" << std::endl;
665
666 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(prec_);
667 if (!prec.is_null()) {
668 prec->setMatrix(tA);
669 reusePreconditioner = true;
670 } else {
671 this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupBlockRelaxation(): reuse of this type is not available (failed cast to CanChangeMatrix), "
672 "reverting to full construction" << std::endl;
673 }
674 }
675
676 if (!reusePreconditioner) {
677 ParameterList& myparamList = const_cast<ParameterList&>(this->GetParameterList());
678 myparamList.print();
679 if(myparamList.isParameter("partitioner: type") &&
680 myparamList.get<std::string>("partitioner: type") == "line") {
681 Teuchos::RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO> > xCoordinates =
683 Teuchos::RCP<Tpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO> > coordinates = Teuchos::rcpFromRef(Xpetra::toTpetra<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO>(*xCoordinates));
684
685 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
686 size_t lclSize = A_->getRangeMap()->getLocalNumElements();
687 if (!matA.is_null())
688 lclSize = matA->getLocalNumRows();
689 size_t numDofsPerNode = lclSize / xCoordinates->getMap()->getLocalNumElements();
690 myparamList.set("partitioner: coordinates", coordinates);
691 myparamList.set("partitioner: PDE equations", (int) numDofsPerNode);
692 }
693
696 prec_->initialize();
697 }
698
699 prec_->compute();
700 }
701
702 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
704 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
705 using STS = Teuchos::ScalarTraits<SC>;
706 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
707 if (!bA.is_null())
708 A_ = bA->Merge();
709
710 RCP<const tRowMatrix> tA = Utilities::Op2NonConstTpetraRow(A_);
711
712 bool reusePreconditioner = false;
713
714 if (this->IsSetup() == true) {
715 // Reuse the constructed preconditioner
716 this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupChebyshev(): Setup() has already been called, assuming reuse" << std::endl;
717
718 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(prec_);
719 if (!prec.is_null()) {
720 prec->setMatrix(tA);
721 reusePreconditioner = true;
722 } else {
723 this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupChebyshev(): reuse of this type is not available (failed cast to CanChangeMatrix), "
724 "reverting to full construction" << std::endl;
725 }
726 }
727
728 // Take care of the eigenvalues
729 ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
730 SC negone = -STS::one();
731 SC lambdaMax = SetupChebyshevEigenvalues(currentLevel,"A","",paramList);
732
733
734 if (!reusePreconditioner) {
737 {
738 SubFactoryMonitor(*this, "Preconditioner init", currentLevel);
739 prec_->initialize();
740 }
741 } else
743
744 {
745 SubFactoryMonitor(*this, "Preconditioner compute", currentLevel);
746 prec_->compute();
747 }
748
749 if (lambdaMax == negone) {
750 typedef Tpetra::RowMatrix<SC, LO, GO, NO> MatrixType;
751
752 Teuchos::RCP<Ifpack2::Chebyshev<MatrixType> > chebyPrec = rcp_dynamic_cast<Ifpack2::Chebyshev<MatrixType> >(prec_);
753 if (chebyPrec != Teuchos::null) {
754 lambdaMax = chebyPrec->getLambdaMaxForApply();
755 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
756 if (!matA.is_null())
757 matA->SetMaxEigenvalueEstimate(lambdaMax);
758 this->GetOStream(Statistics1) << "chebyshev: max eigenvalue (calculated by Ifpack2)" << " = " << lambdaMax << std::endl;
759 }
760 TEUCHOS_TEST_FOR_EXCEPTION(lambdaMax == negone, Exceptions::RuntimeError, "MueLu::Ifpack2Smoother::Setup(): no maximum eigenvalue estimate");
761 }
762 }
763
764
765 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
767 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
768 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
769 if (!bA.is_null())
770 A_ = bA->Merge();
771
772 RCP<const tRowMatrix> tA = Utilities::Op2NonConstTpetraRow(A_);
773
774 bool reusePreconditioner = false;
775 if (this->IsSetup() == true) {
776 // Reuse the constructed preconditioner
777 this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupHiptmair(): Setup() has already been called, assuming reuse" << std::endl;
778
779 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(prec_);
780 if (!prec.is_null()) {
781 prec->setMatrix(tA);
782 reusePreconditioner = true;
783 } else {
784 this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupHiptmair(): reuse of this type is not available (failed cast to CanChangeMatrix), "
785 "reverting to full construction" << std::endl;
786 }
787 }
788
789 // If we're doing Chebyshev subsmoothing, we'll need to get the eigenvalues
790 ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
791 std::string smoother1 = paramList.get("hiptmair: smoother type 1","CHEBYSHEV");
792 std::string smoother2 = paramList.get("hiptmair: smoother type 2","CHEBYSHEV");
793 // SC lambdaMax11,lambdaMax22;
794
795 if(smoother1 == "CHEBYSHEV") {
796 ParameterList & list1 = paramList.sublist("hiptmair: smoother list 1");
797 //lambdaMax11 =
798 SetupChebyshevEigenvalues(currentLevel,"A","EdgeMatrix ",list1);
799 }
800 if(smoother2 == "CHEBYSHEV") {
801 ParameterList & list2 = paramList.sublist("hiptmair: smoother list 2");
802 //lambdaMax22 =
803 SetupChebyshevEigenvalues(currentLevel,"A","EdgeMatrix ",list2);
804 }
805
806 // FIXME: Should really add some checks to make sure the eigenvalue calcs worked like in
807 // the regular SetupChebyshev
808
809 // Grab the auxillary matrices and stick them on the list
810 RCP<Operator> NodeMatrix = currentLevel.Get<RCP<Operator> >("NodeMatrix");
811 RCP<Operator> D0 = currentLevel.Get<RCP<Operator> >("D0");
812
813 RCP<tRowMatrix> tNodeMatrix = Utilities::Op2NonConstTpetraRow(NodeMatrix);
814 RCP<tRowMatrix> tD0 = Utilities::Op2NonConstTpetraRow(D0);
815
816
817 Teuchos::ParameterList newlist;
818 newlist.set("P",tD0);
819 newlist.set("PtAP",tNodeMatrix);
820 if (!reusePreconditioner) {
822 SetPrecParameters(newlist);
823 prec_->initialize();
824 }
825
826 prec_->compute();
827 }
828
829
830
831 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
832 Scalar Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::SetupChebyshevEigenvalues(Level & currentLevel, const std::string & matrixName, const std::string & label, ParameterList & paramList) const {
833 // Helper: This gets used for smoothers that want to set up Chebyhev
834 typedef Teuchos::ScalarTraits<SC> STS;
835 SC negone = -STS::one();
836 RCP<const Matrix> currentA = currentLevel.Get<RCP<Matrix> >(matrixName);
837 SC lambdaMax = negone;
838
839 std::string maxEigString = "chebyshev: max eigenvalue";
840 std::string eigRatioString = "chebyshev: ratio eigenvalue";
841
842 // Get/calculate the maximum eigenvalue
843 if (paramList.isParameter(maxEigString)) {
844 if (paramList.isType<double>(maxEigString))
845 lambdaMax = paramList.get<double>(maxEigString);
846 else
847 lambdaMax = paramList.get<SC>(maxEigString);
848 this->GetOStream(Statistics1) << label << maxEigString << " (cached with smoother parameter list) = " << lambdaMax << std::endl;
849 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
850 if (!matA.is_null())
851 matA->SetMaxEigenvalueEstimate(lambdaMax);
852
853 } else {
854 lambdaMax = currentA->GetMaxEigenvalueEstimate();
855 if (lambdaMax != negone) {
856 this->GetOStream(Statistics1) << label << maxEigString << " (cached with matrix) = " << lambdaMax << std::endl;
857 paramList.set(maxEigString, lambdaMax);
858 }
859 }
860
861 // Calculate the eigenvalue ratio
862 const SC defaultEigRatio = 20;
863
864 SC ratio = defaultEigRatio;
865 if (paramList.isParameter(eigRatioString)) {
866 if (paramList.isType<double>(eigRatioString))
867 ratio = paramList.get<double>(eigRatioString);
868 else
869 ratio = paramList.get<SC>(eigRatioString);
870 }
871 if (currentLevel.GetLevelID()) {
872 // Update ratio to be
873 // ratio = max(number of fine DOFs / number of coarse DOFs, defaultValue)
874 //
875 // NOTE: We don't need to request previous level matrix as we know for sure it was constructed
876 RCP<const Matrix> fineA = currentLevel.GetPreviousLevel()->Get<RCP<Matrix> >(matrixName);
877 size_t nRowsFine = fineA->getGlobalNumRows();
878 size_t nRowsCoarse = currentA->getGlobalNumRows();
879
880 SC levelRatio = as<SC>(as<float>(nRowsFine)/nRowsCoarse);
881 if (STS::magnitude(levelRatio) > STS::magnitude(ratio))
882 ratio = levelRatio;
883 }
884
885 this->GetOStream(Statistics1) << label << eigRatioString << " (computed) = " << ratio << std::endl;
886 paramList.set(eigRatioString, ratio);
887
888 if (paramList.isParameter("chebyshev: use rowsumabs diagonal scaling")) {
889 this->GetOStream(Runtime1) << "chebyshev: using rowsumabs diagonal scaling" << std::endl;
890 bool doScale = false;
891 doScale = paramList.get<bool>("chebyshev: use rowsumabs diagonal scaling");
892 paramList.remove("chebyshev: use rowsumabs diagonal scaling");
893 double chebyReplaceTol = Teuchos::ScalarTraits<Scalar>::eps()*100;
894 std::string paramName = "chebyshev: rowsumabs diagonal replacement tolerance";
895 if (paramList.isParameter(paramName)) {
896 chebyReplaceTol = paramList.get<double>(paramName);
897 paramList.remove(paramName);
898 }
899 double chebyReplaceVal = Teuchos::ScalarTraits<double>::zero();
900 paramName = "chebyshev: rowsumabs diagonal replacement value";
901 if (paramList.isParameter(paramName)) {
902 chebyReplaceVal = paramList.get<double>(paramName);
903 paramList.remove(paramName);
904 }
905 bool chebyReplaceSingleEntryRowWithZero = false;
906 paramName = "chebyshev: rowsumabs replace single entry row with zero";
907 if (paramList.isParameter(paramName)) {
908 chebyReplaceSingleEntryRowWithZero = paramList.get<bool>(paramName);
909 paramList.remove(paramName);
910 }
911 bool useAverageAbsDiagVal = false;
912 paramName = "chebyshev: rowsumabs use automatic diagonal tolerance";
913 if (paramList.isParameter(paramName)) {
914 useAverageAbsDiagVal = paramList.get<bool>(paramName);
915 paramList.remove(paramName);
916 }
917 if (doScale) {
918 const bool doReciprocal = true;
919 RCP<Vector> lumpedDiagonal = Utilities::GetLumpedMatrixDiagonal(*currentA, doReciprocal, chebyReplaceTol, chebyReplaceVal, chebyReplaceSingleEntryRowWithZero, useAverageAbsDiagVal);
920 const Xpetra::TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmpVec = dynamic_cast<const Xpetra::TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>&>(*lumpedDiagonal);
921 paramList.set("chebyshev: operator inv diagonal",tmpVec.getTpetra_Vector());
922 }
923 }
924
925 return lambdaMax;
926 }
927
928
929
930
931 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
933 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
934 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
935 if (!bA.is_null())
936 A_ = bA->Merge();
937
938 RCP<const tRowMatrix> tA = Utilities::Op2NonConstTpetraRow(A_);
939
940 bool reusePreconditioner = false;
941 if (this->IsSetup() == true) {
942 // Reuse the constructed preconditioner
943 this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupGeneric(): Setup() has already been called, assuming reuse" << std::endl;
944
945 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(prec_);
946 if (!prec.is_null()) {
947 prec->setMatrix(tA);
948 reusePreconditioner = true;
949 } else {
950 this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupGeneric(): reuse of this type is not available (failed cast to CanChangeMatrix), "
951 "reverting to full construction" << std::endl;
952 }
953 }
954
955 if (!reusePreconditioner) {
958 prec_->initialize();
959 }
960
961 prec_->compute();
962 }
963
964 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
965 void Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector& X, const MultiVector& B, bool InitialGuessIsZero) const {
966 TEUCHOS_TEST_FOR_EXCEPTION(SmootherPrototype::IsSetup() == false, Exceptions::RuntimeError, "MueLu::Ifpack2Smoother::Apply(): Setup() has not been called");
967
968 // Forward the InitialGuessIsZero option to Ifpack2
969 // TODO: It might be nice to switch back the internal
970 // "zero starting solution" option of the ifpack2 object prec_ to his
971 // initial value at the end but there is no way right now to get
972 // the current value of the "zero starting solution" in ifpack2.
973 // It's not really an issue, as prec_ can only be used by this method.
974 Teuchos::ParameterList paramList;
975 bool supportInitialGuess = false;
976 const Teuchos::ParameterList params = this->GetParameterList();
977
978 if (prec_->supportsZeroStartingSolution()) {
979 prec_->setZeroStartingSolution(InitialGuessIsZero);
980 supportInitialGuess = true;
981 } else if (type_ == "SCHWARZ") {
982 paramList.set("schwarz: zero starting solution", InitialGuessIsZero);
983 //Because additive Schwarz has "delta" semantics, it's sufficient to
984 //toggle only the zero initial guess flag, and not pass in already
985 //set parameters. If we call SetPrecParameters, the subdomain solver
986 //will be destroyed.
987 prec_->setParameters(paramList);
988 supportInitialGuess = true;
989 }
990
991 //TODO JJH 30Apr2014 Calling SetPrecParameters(paramList) when the smoother
992 //is Ifpack2::AdditiveSchwarz::setParameterList() will destroy the subdomain
993 //(aka inner) solver. This behavior is documented but a departure from what
994 //it previously did, and what other Ifpack2 solvers currently do. So I have
995 //moved SetPrecParameters(paramList) into the if-else block above.
996
997 // Apply
998 if (InitialGuessIsZero || supportInitialGuess) {
999 Tpetra::MultiVector<SC,LO,GO,NO>& tpX = Utilities::MV2NonConstTpetraMV(X);
1000 const Tpetra::MultiVector<SC,LO,GO,NO>& tpB = Utilities::MV2TpetraMV(B);
1001 prec_->apply(tpB, tpX);
1002 } else {
1003 typedef Teuchos::ScalarTraits<Scalar> TST;
1004
1005 RCP<MultiVector> Residual;
1006 {
1007 std::string prefix = this->ShortClassName() + ": Apply: ";
1008 RCP<TimeMonitor> tM = rcp(new TimeMonitor(*this, prefix + "residual calculation", Timings0));
1009 Residual = Utilities::Residual(*A_, X, B);
1010 }
1011
1012 RCP<MultiVector> Correction = MultiVectorFactory::Build(A_->getDomainMap(), X.getNumVectors());
1013
1014 Tpetra::MultiVector<SC,LO,GO,NO>& tpX = Utilities::MV2NonConstTpetraMV(*Correction);
1015 const Tpetra::MultiVector<SC,LO,GO,NO>& tpB = Utilities::MV2TpetraMV(*Residual);
1016
1017 prec_->apply(tpB, tpX);
1018
1019 X.update(TST::one(), *Correction, TST::one());
1020 }
1021 }
1022
1023 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
1024 RCP<MueLu::SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Copy() const {
1025 RCP<Ifpack2Smoother> smoother = rcp(new Ifpack2Smoother(*this) );
1026 smoother->SetParameterList(this->GetParameterList());
1027 return smoother;
1028 }
1029
1030 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
1032 std::ostringstream out;
1034 out << prec_->description();
1035 } else {
1037 out << "{type = " << type_ << "}";
1038 }
1039 return out.str();
1040 }
1041
1042 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
1043 void Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::print(Teuchos::FancyOStream &out, const VerbLevel verbLevel) const {
1045
1046 if (verbLevel & Parameters0)
1047 out0 << "Prec. type: " << type_ << std::endl;
1048
1049 if (verbLevel & Parameters1) {
1050 out0 << "Parameter list: " << std::endl;
1051 Teuchos::OSTab tab2(out);
1052 out << this->GetParameterList();
1053 out0 << "Overlap: " << overlap_ << std::endl;
1054 }
1055
1056 if (verbLevel & External)
1057 if (prec_ != Teuchos::null) {
1058 Teuchos::OSTab tab2(out);
1059 out << *prec_ << std::endl;
1060 }
1061
1062 if (verbLevel & Debug) {
1063 out0 << "IsSetup: " << Teuchos::toString(SmootherPrototype::IsSetup()) << std::endl
1064 << "-" << std::endl
1065 << "RCP<prec_>: " << prec_ << std::endl;
1066 }
1067 }
1068
1069 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
1071 typedef Tpetra::RowMatrix<SC,LO,GO,NO> MatrixType;
1072 // NOTE: Only works for a subset of Ifpack2's smoothers
1073 RCP<Ifpack2::Relaxation<MatrixType> > pr = rcp_dynamic_cast<Ifpack2::Relaxation<MatrixType> >(prec_);
1074 if(!pr.is_null()) return pr->getNodeSmootherComplexity();
1075
1076 RCP<Ifpack2::Chebyshev<MatrixType> > pc = rcp_dynamic_cast<Ifpack2::Chebyshev<MatrixType> >(prec_);
1077 if(!pc.is_null()) return pc->getNodeSmootherComplexity();
1078
1079 RCP<Ifpack2::BlockRelaxation<MatrixType> > pb = rcp_dynamic_cast<Ifpack2::BlockRelaxation<MatrixType> >(prec_);
1080 if(!pb.is_null()) return pb->getNodeSmootherComplexity();
1081
1082 RCP<Ifpack2::ILUT<MatrixType> > pi = rcp_dynamic_cast<Ifpack2::ILUT<MatrixType> >(prec_);
1083 if(!pi.is_null()) return pi->getNodeSmootherComplexity();
1084
1085 RCP<Ifpack2::RILUK<MatrixType> > pk = rcp_dynamic_cast<Ifpack2::RILUK<MatrixType> >(prec_);
1086 if(!pk.is_null()) return pk->getNodeSmootherComplexity();
1087
1088
1089 return Teuchos::OrdinalTraits<size_t>::invalid();
1090 }
1091
1092
1093} // namespace MueLu
1094
1095#endif // HAVE_MUELU_IFPACK2
1096#endif // MUELU_IFPACK2SMOOTHER_DEF_HPP
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
static Teuchos::RCP< Preconditioner< typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > create(const std::string &precType, const Teuchos::RCP< const MatrixType > &matrix)
virtual std::string ShortClassName() const
Return the class name of the object, without template parameters and without namespace.
virtual std::string description() const
Return a simple one-line description of this object.
Exception indicating invalid cast attempted.
Exception throws to report errors in the internal logical of the program.
Timer to be used in factories. Similar to Monitor but with additional timers.
void Input(Level &level, const std::string &varName) const
T Get(Level &level, const std::string &varName) const
RCP< ParameterList > RemoveFactoriesFromList(const ParameterList &list) const
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
void SetupSchwarz(Level &currentLevel)
void Setup(Level &currentLevel)
Set up the smoother.
Scalar SetupChebyshevEigenvalues(Level &currentLevel, const std::string &matrixName, const std::string &label, Teuchos::ParameterList &paramList) const
void SetupAggregate(Level &currentLevel)
void SetupBlockRelaxation(Level &currentLevel)
void SetupChebyshev(Level &currentLevel)
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
void SetParameterList(const Teuchos::ParameterList &paramList)
Set parameters from a parameter list and return with default values.
RCP< SmootherPrototype > Copy() const
std::string description() const
Return a simple one-line description of this object.
void SetupHiptmair(Level &currentLevel)
LO overlap_
overlap when using the smoother in additive Schwarz mode
void SetPrecParameters(const Teuchos::ParameterList &list=Teuchos::ParameterList()) const
void SetupTopological(Level &currentLevel)
RCP< Ifpack2::Preconditioner< Scalar, LocalOrdinal, GlobalOrdinal, Node > > prec_
pointer to Ifpack2 preconditioner object
std::string type_
ifpack2-specific key phrase that denote smoother type
friend class Ifpack2Smoother
Constructor.
void DeclareInput(Level &currentLevel) const
Input.
RCP< Operator > A_
matrix, used in apply if solving residual equation
void SetupLineSmoothing(Level &currentLevel)
void SetupGeneric(Level &currentLevel)
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the preconditioner.
Class that holds all level-specific information.
RCP< Level > & GetPreviousLevel()
Previous level.
int GetLevelID() const
Return level number.
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access)....
virtual const Teuchos::ParameterList & GetParameterList() const
virtual void SetParameterList(const Teuchos::ParameterList &paramList)=0
Set parameters from a parameter list and return with default values.
void declareConstructionOutcome(bool fail, std::string msg)
bool IsSetup() const
Get the state of a smoother prototype.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
Integrates Teuchos::TimeMonitor with MueLu verbosity system.
static Teuchos::ArrayRCP< const bool > DetectDirichletRows(const Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &A, const Magnitude &tol=Teuchos::ScalarTraits< Magnitude >::zero(), bool count_twos_as_dirichlet=false)
static Teuchos::RCP< Vector > GetLumpedMatrixDiagonal(Matrix const &A, const bool doReciprocal=false, Magnitude tol=Teuchos::ScalarTraits< Scalar >::magnitude(Teuchos::ScalarTraits< Scalar >::zero()), Scalar valReplacement=Teuchos::ScalarTraits< Scalar >::zero(), const bool replaceSingleEntryRowWithZero=false, const bool useAverageAbsDiagVal=false)
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &Op, const MultiVector &X, const MultiVector &RHS)
static RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraRow(RCP< Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2TpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > const vec)
Helper utility to pull out the underlying Tpetra objects from an Xpetra object.
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec)
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
bool IsPrint(MsgType type, int thisProcRankOnly=-1) const
Find out whether we need to print out information for a specific message type.
void FindGeometricSeedOrdinals(Teuchos::RCP< Basis > basis, const LOFieldContainer &elementToNodeMap, std::vector< std::vector< LocalOrdinal > > &seeds, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &rowMap, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &columnMap)
Teuchos::RCP< Intrepid2::Basis< KokkosExecutionSpace, Scalar, Scalar > > BasisFactory(const std::string &name, int &degree)
Namespace for MueLu classes and methods.
@ Warnings0
Important warning messages (one line).
@ Debug
Print additional debugging information.
@ Statistics1
Print more statistics.
@ External
Print external lib objects.
@ Timings0
High level timing information (use Teuchos::TimeMonitor::summarize() to print).
@ Runtime0
One-liner description of what is happening.
@ Runtime1
Description of what is happening (more verbose).
@ Timings
Print all timing information.
@ Parameters0
Print class parameters.
@ Statistics0
Print statistics that do not involve significant additional computation.
@ Parameters1
Print class parameters (more parameters, more verbose).