LifeV
MultiscaleAlgorithmAitken.cpp
Go to the documentation of this file.
1 //@HEADER
2 /*
3 *******************************************************************************
4 
5  Copyright (C) 2004, 2005, 2007 EPFL, Politecnico di Milano, INRIA
6  Copyright (C) 2010 EPFL, Politecnico di Milano, Emory University
7 
8  This file is part of LifeV.
9 
10  LifeV is free software; you can redistribute it and/or modify
11  it under the terms of the GNU Lesser General Public License as published by
12  the Free Software Foundation, either version 3 of the License, or
13  (at your option) any later version.
14 
15  LifeV is distributed in the hope that it will be useful,
16  but WITHOUT ANY WARRANTY; without even the implied warranty of
17  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18  Lesser General Public License for more details.
19 
20  You should have received a copy of the GNU Lesser General Public License
21  along with LifeV. If not, see <http://www.gnu.org/licenses/>.
22 
23 *******************************************************************************
24 */
25 //@HEADER
26 
27 /*!
28  * @file
29  * @brief File containing the Multiscale Aitken Algorithm
30  *
31  * @date 23-10-2009
32  * @author Cristiano Malossi <cristiano.malossi@epfl.ch>
33  *
34  * @maintainer Cristiano Malossi <cristiano.malossi@epfl.ch>
35  */
36 
37 #include <lifev/multiscale/algorithms/MultiscaleAlgorithmAitken.hpp>
38 
39 namespace LifeV
40 {
41 namespace Multiscale
42 {
43 
44 // ===================================================
45 // Constructors & Destructor
46 // ===================================================
49  M_methodMap (),
50  M_method (),
52 {
53 
54 #ifdef HAVE_LIFEV_DEBUG
55  debugStream ( 8011 ) << "MultiscaleAlgorithmAitken::MultiscaleAlgorithmAitken() \n";
56 #endif
57 
58  M_type = Aitken;
59 
60  //Set Map
61  M_methodMap["Scalar"] = Scalar;
62  M_methodMap["Vectorial"] = Vectorial;
63  //M_methodMap["VectorialBlock"] = VectorialBlock;
64 }
65 
66 // ===================================================
67 // Multiscale Algorithm Virtual Methods
68 // ===================================================
69 void
70 MultiscaleAlgorithmAitken::setupData ( const std::string& fileName )
71 {
72 
73 #ifdef HAVE_LIFEV_DEBUG
74  debugStream ( 8011 ) << "MultiscaleAlgorithmAitken::setupData( fileName ) \n";
75 #endif
76 
77  // Read parameters
78  multiscaleParameterListPtr_Type solverParametersList = Teuchos::rcp ( new Teuchos::ParameterList );
79  solverParametersList = Teuchos::getParametersFromXmlFile ( fileName );
80 
81  // Set main parameters
82  setAlgorithmName ( solverParametersList->sublist ( "Multiscale", true, "" ) );
83  setAlgorithmParameters ( solverParametersList->sublist ( "Multiscale Algorithm", true, "" ) );
84 }
85 
86 void
87 MultiscaleAlgorithmAitken::subIterate()
88 {
89 
90 #ifdef HAVE_LIFEV_DEBUG
91  debugStream ( 8011 ) << "MultiscaleAlgorithmAitken::subIterate() \n";
92 #endif
93 
94  multiscaleAlgorithm_Type::subIterate();
95 
96  // Verify tolerance
97  if ( checkResidual ( 0 ) )
98  {
99  return;
100  }
101 
102  M_multiscale->exportCouplingVariables ( *M_couplingVariables );
103 
104  M_generalizedAitken.restart();
105 
106  // Temporary Computation of a Block Vector - Testing purpose
107  // VectorType blocksVector( M_couplingVariables ); blocksVector = 0.0;
108  // for ( UInt i = 1 ; i < blocksVector.size() ; i = i+2)
109  // blocksVector[i] = 1.0;
110  // std::cout << "blocksVector: " << std::endl;
111  // blocksVector.showMe();
112 
113  for ( UInt subIT = 1; subIT <= M_subiterationsMaximumNumber; ++subIT )
114  {
115  // Update Coupling Variables
116  switch ( M_method )
117  {
118  case Scalar:
119 
120  *M_couplingVariables += M_generalizedAitken.computeDeltaLambdaScalar ( *M_couplingVariables, *M_couplingResiduals );
121 
122  break;
123 
124  case Vectorial:
125 
126  *M_couplingVariables += M_generalizedAitken.computeDeltaLambdaVector ( *M_couplingVariables, *M_couplingResiduals, true );
127 
128  break;
129 
130  case VectorialBlock:
131 
132  //*M_couplingVariables += M_generalizedAitken.computeDeltaLambdaVectorBlock( *M_couplingVariables, *M_couplingResiduals, blocksVector, 2 );
133 
134  break;
135  }
136 
137  // Import Coupling Variables inside the coupling blocks
138  M_multiscale->importCouplingVariables ( *M_couplingVariables );
139 
140  // Verify tolerance
141  if ( checkResidual ( subIT ) )
142  {
143  return;
144  }
145  }
146 
147  save ( M_subiterationsMaximumNumber, M_couplingResiduals->norm2() );
148 
149  multiscaleErrorCheck ( Tolerance, "Aitken algorithm residual: " + number2string ( M_couplingResiduals->norm2() ) +
150  " (required: " + number2string ( M_tolerance ) + ")\n", M_multiscale->communicator() == 0 );
151 }
152 
153 void
154 MultiscaleAlgorithmAitken::showMe()
155 {
156  if ( M_comm->MyPID() == 0 )
157  {
158  multiscaleAlgorithm_Type::showMe();
159 
160  std::cout << "Aitken Method = " << enum2String ( M_method, M_methodMap ) << std::endl;
161  std::cout << std::endl << std::endl;
162  }
163 }
164 
165 // ===================================================
166 // Set Methods
167 // ===================================================
168 void
169 MultiscaleAlgorithmAitken::setAlgorithmParameters ( const multiscaleParameterList_Type& parameterList )
170 {
171  multiscaleAlgorithm_Type::setAlgorithmParameters ( parameterList );
172 
173  M_generalizedAitken.setDefaultOmega ( parameterList.get<Real> ( "Omega" ) );
174  M_generalizedAitken.useDefaultOmega ( parameterList.get<bool> ( "Fixed omega" ) );
175  M_generalizedAitken.setOmegaMin ( parameterList.get<Real> ( "Range minimum" ) );
176  M_generalizedAitken.setOmegaMax ( parameterList.get<Real> ( "Range maximum" ) );
177  M_generalizedAitken.setMinimizationType ( parameterList.get<bool> ( "Inverse omega" ) );
178  M_method = M_methodMap[ parameterList.get<std::string> ( "Method" ) ];
179 }
180 
181 } // Namespace Multiscale
182 } // Namespace LifeV
MultiscaleAlgorithmAitken - The Multiscale Algorithm implementation of Aitken.
void setupData(const std::string &fileName)
Setup the data of the algorithm using a data file.
MultiscaleAlgorithm multiscaleAlgorithm_Type
LinearSolver::parameterListPtr_Type multiscaleParameterListPtr_Type