Nektar++
DriverArpack.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File DriverArpack.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Arnoldi solver using Arpack
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 
38 namespace Nektar
39 {
40  namespace SolverUtils
41  {
43  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","LargestReal" ,0),
44  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","SmallestReal" ,1),
45  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","LargestImag" ,2),
46  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","SmallestImag" ,3),
47  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","LargestMag" ,4),
48  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","SmallestMag" ,5),
49  };
52 
54 
56  { "LR", "SR", "LI", "SI", "LM", "SM" };
57 
58 
59  /**
60  *
61  */
63  : DriverArnoldi(pSession)
64  {
65  }
66 
67 
68  /**
69  *
70  */
72  {
73  }
74 
75  /**
76  *
77  */
78  void DriverArpack::v_InitObject(ostream &out)
79  {
81 
82  //Initialisation of Arnoldi parameters
83  m_maxn = 1000000; // Maximum size of the problem
84  m_maxnev = 200; // maximum number of eigenvalues requested
85  m_maxncv = 500; // Largest number of basis vector used in Implicitly Restarted Arnoldi
86 
87  // Error alerts
88  ASSERTL0(m_nvec <= m_maxnev,"NEV is greater than MAXNEV");
89  ASSERTL0(m_kdim <= m_maxncv,"NEV is greater than MAXNEV");
90  ASSERTL0(2 <= m_kdim-m_nvec,"NCV-NEV is less than 2");
91 
92  m_equ[0]->PrintSummary(out);
93 
94  // Print session parameters
95  out << "\tArnoldi solver type : Arpack" << endl;
96 
97  out << "\tArpack problem type : ";
98  out << ArpackProblemTypeTrans[m_session->GetSolverInfoAsEnum<int>("ArpackProblemType")] << endl;
100 
101  m_equ[m_nequ - 1]->DoInitialise();
102 
103  //FwdTrans Initial conditions to be in Coefficient Space
104  m_equ[m_nequ-1] ->TransPhysToCoeff();
105 
106  }
107 
108 
109  void DriverArpack::v_Execute(ostream &out)
110 
111  {
112  Array<OneD, NekDouble> tmpworkd;
113 
114  int nq = m_equ[0]->UpdateFields()[0]->GetNcoeffs(); // Number of points in the mesh
115  int n = m_nfields*nq; // Number of points in eigenvalue calculation
116  int lworkl = 3*m_kdim*(m_kdim+2); // Size of work array
117  int ido ; //REVERSE COMMUNICATION parameter. At the first call must be initialised at 0
118  int info; // do not set initial vector (info=0 random initial vector, info=1 read initial vector from session file)
119 
120  int iparam[11];
121  int ipntr[14];
122 
123  Array<OneD, int> ritzSelect;
126  Array<OneD, NekDouble> workev;
128  NekDouble sigmar, sigmai;
129 
130  Array<OneD, NekDouble> resid(n);
132  Array<OneD, NekDouble> workl(lworkl);
133  Array<OneD, NekDouble> workd(3*n, 0.0);
134 
135  ASSERTL0(n <= m_maxn, "N is greater than MAXN");
136 
137  if(m_session->DefinesFunction("InitialConditions"))
138  {
139 
140  out << "\tInital vector : input file " << endl;
141  info = 1;
143 
144  }
145  else
146  {
147  out << "\tInital vector : random " << endl;
148  info = 0;
149  }
150 
151 
152  iparam[0] = 1; // strategy for shift-invert
153  iparam[1] = 0; // (deprecated)
154  iparam[2] = m_nits; // maximum number of iterations allowed/taken
155  iparam[3] = 1; // blocksize to be used for recurrence
156  iparam[4] = 0; // number of converged ritz eigenvalues
157  iparam[5] = 0; // (deprecated)
158  if((fabs(m_realShift) > NekConstants::kNekZeroTol)|| // use shift if m_realShift > 1e-12
160  {
161  iparam[6] = 3;
162  }
163  else
164  {
165  iparam[6] = 1; // computation mode 1=> matrix-vector prod
166  }
167  iparam[7] = 0; // (for shift-invert)
168  iparam[8] = 0; // number of MV operations
169  iparam[9] = 0; // number of BV operations
170  iparam[10]= 0; // number of reorthogonalisation steps
171 
172  int cycle = 0;
173  const char* problem = ArpackProblemTypeTrans[m_session->GetSolverInfoAsEnum<int>("ArpackProblemType")].c_str();
174 
175  std::string name = m_session->GetSessionName() + ".evl";
176  ofstream pFile(name.c_str());
177 
178  ido = 0; //At the first call must be initialisedat 0
179 
180  while(ido != 99)//ido==-1 || ido==1 || ido==0)
181  {
182  //Routine for eigenvalue evaluation for non-symmetric operators
183  Arpack::Dnaupd( ido, "I", // B='I' for std eval problem
184  n, problem, m_nvec,
185  m_evtol, &resid[0], m_kdim,
186  &v[0], n, iparam, ipntr, &workd[0],
187  &workl[0], lworkl, info);
188 
189  //Plotting of real and imaginary part of the eigenvalues from workl
190  out << "\rIteration " << cycle << ", output: " << info << ", ido=" << ido << " " << std::flush;
191 
192  if(!((cycle-1)%m_kdim)&&(cycle> m_kdim))
193  {
194  pFile << "Krylov spectrum at iteration: " << cycle << endl;
195 
197  {
198  pFile << "EV Magnitude Angle Growth Frequency Residual" << endl;
199  }
200  else
201  {
202  pFile << "EV Real Imaginary inverse real inverse imag Residual" << endl;
203  }
204 
205  out << endl;
206  for(int k=0; k<=m_kdim-1; ++k)
207  {
208  // write m_nvec eigs to screen
209  if(m_kdim-1-k < m_nvec)
210  {
211  WriteEvs(out,k, workl[ipntr[5]-1+k],workl[ipntr[6]-1+k]);
212  }
213  // write m_kdim eigs to screen
214  WriteEvs(pFile,k, workl[ipntr[5]-1+k],workl[ipntr[6]-1+k]);
215  }
216  }
217 
218  cycle++;
219 
220  if (ido == 99) break;
221 
222  ASSERTL0(ido == 1, "Unexpected reverse communication request.");
223 
224  //workd[inptr[0]-1] copied into operator fields
225  CopyArnoldiArrayToField(tmpworkd = workd + (ipntr[0]-1));
226 
227  m_equ[0]->TransCoeffToPhys();
228 
229  m_equ[0]->DoSolve();
230 
231  if(!(cycle%m_infosteps))
232  {
233  out << endl;
234  m_equ[0]->Output();
235  }
236 
238  {
239  //start Adjoint with latest fields of direct
240  CopyFwdToAdj();
241 
242  m_equ[1]->TransCoeffToPhys();
243  m_equ[1]->DoSolve();
244  }
245 
246  // operated fields are copied into workd[inptr[1]-1]
247  CopyFieldToArnoldiArray(tmpworkd = workd + (ipntr[1]-1));
248 
249  }
250 
251  out<< endl << "Converged in " << iparam[8] << " iterations" << endl;
252 
253  ASSERTL0(info >= 0," Error with Dnaupd");
254 
255  ritzSelect = Array<OneD, int> (m_kdim,0);
256  dr = Array<OneD, NekDouble> (m_nvec+1,0.0);
257  di = Array<OneD, NekDouble> (m_nvec+1,0.0);
258  workev = Array<OneD, NekDouble> (3*m_kdim);
259  z = Array<OneD, NekDouble> (n*(m_nvec+1));
260 
261  sigmar = m_realShift;
262  sigmai = m_imagShift;
263 
264  //Setting 'A', Ritz vectors are computed. 'S' for Shur vectors
265  Arpack::Dneupd(1, "A", ritzSelect.get(), dr.get(), di.get(), z.get(), n, sigmar, sigmai, workev.get(), "I", n, problem, m_nvec, m_evtol, resid.get(), m_kdim, v.get(), n, iparam, ipntr, workd.get(), workl.get(),lworkl,info);
266 
267  ASSERTL0(info == 0, " Error with Dneupd");
268  int nconv=iparam[4];
269  Array<OneD, MultiRegions::ExpListSharedPtr> fields = m_equ[0]->UpdateFields();
270 
271  out << "Converged Eigenvalues: " << nconv << endl;
272  pFile << "Converged Eigenvalues:"<< nconv << endl;
273 
275  {
276  pFile << "EV Magnitude Angle Growth Frequency" << endl;
277  }
278  else
279  {
280  pFile << "EV Real Imaginary inverse real inverse imag" << endl;
281  }
282 
283 
284  for(int i= 0; i< nconv; ++i)
285  {
286  WriteEvs(out,i,dr[i],di[i]);
287  WriteEvs(pFile,i,dr[i],di[i]);
288 
289  std::string file = m_session->GetSessionName() + "_eig_"
290  + boost::lexical_cast<std::string>(i);
291  WriteFld(file,z + i*nq);
292  }
293 
294  m_real_evl = dr;
295  m_imag_evl = di;
296 
297  pFile.close();
298 
299  for(int j = 0; j < m_equ[0]->GetNvariables(); ++j)
300  {
301  NekDouble vL2Error = m_equ[0]->L2Error(j,false);
302  NekDouble vLinfError = m_equ[0]->LinfError(j);
303  if (m_comm->GetRank() == 0)
304  {
305  out << "L 2 error (variable " << m_equ[0]->GetVariable(j) << ") : " << vL2Error << endl;
306  out << "L inf error (variable " << m_equ[0]->GetVariable(j) << ") : " << vLinfError << endl;
307  }
308  }
309  }
310  }
311 }
Array< OneD, NekDouble > m_imag_evl
Definition: DriverArnoldi.h:68
static void Dneupd(const int &rvec, const char *howmny, const int *select, double *dr, double *di, double *z, const int &ldz, const double &sigmar, const double &sigmai, double *workev, const char *bmat, const int &n, const char *which, const int &nev, const double &tol, double *resid, const int &ncv, double *v, const int &ldv, int *iparam, int *ipntr, double *workd, double *workl, const int &lworkl, int &info)
Post-processing routine to computed eigenvector of computed eigenvalues in Dnaupd.
Definition: Arpack.hpp:149
virtual ~DriverArpack()
Destructor.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
void CopyFwdToAdj()
Copy the forward field to the adjoint system in transient growth calculations.
NekDouble m_evtol
Maxmum number of iterations.
Definition: DriverArnoldi.h:57
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
void CopyArnoldiArrayToField(Array< OneD, NekDouble > &array)
Copy Arnoldi storage to fields.
SOLVER_UTILS_EXPORT void ArnoldiSummary(std::ostream &out)
void CopyFieldToArnoldiArray(Array< OneD, NekDouble > &array)
Copy fields to Arnoldi storage.
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:50
void WriteFld(std::string file, std::vector< Array< OneD, NekDouble > > coeffs)
static std::string arpackProblemTypeLookupIds[]
Definition: DriverArpack.h:84
void WriteEvs(ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble)
DriverArpack(const LibUtilities::SessionReaderSharedPtr pSession)
Constructor.
static const NekDouble kNekZeroTol
int m_nvec
Dimension of Krylov subspace.
Definition: DriverArnoldi.h:55
Array< OneD, NekDouble > m_real_evl
Definition: DriverArnoldi.h:67
static void Dnaupd(int &ido, const char *bmat, const int &n, const char *which, const int &nev, const double &tol, double *resid, const int &ncv, double *v, const int &ldv, int *iparam, int *ipntr, double *workd, double *workl, const int &lworkl, int &info)
Top level reverse communication interface to solve real double-precision non-symmetric problems...
Definition: Arpack.hpp:135
enum EvolutionOperatorType m_EvolutionOperator
Evolution Operator.
Definition: Driver.h:100
double NekDouble
int m_nfields
interval to dump information if required.
Definition: DriverArnoldi.h:63
bool m_timeSteppingAlgorithm
Period of time stepping algorithm.
Definition: DriverArnoldi.h:59
static std::string ArpackProblemTypeTrans[]
Definition: DriverArpack.h:89
static std::string driverLookupId
Definition: DriverArpack.h:86
static std::string arpackProblemTypeDefault
Definition: DriverArpack.h:85
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:94
int m_infosteps
underlying operator is time stepping
Definition: DriverArnoldi.h:61
DriverFactory & GetDriverFactory()
Definition: Driver.cpp:64
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:85
static std::string RegisterDefaultSolverInfo(const std::string &pName, const std::string &pValue)
Registers the default string value of a solver info property.
static DriverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
Creates an instance of this class.
Definition: DriverArpack.h:54
virtual void v_Execute(ostream &out=cout)
Virtual function for solve implementation.
virtual void v_InitObject(ostream &out=cout)
Virtual function for initialisation implementation.
Base class for the development of solvers.
Definition: DriverArnoldi.h:46
virtual void v_InitObject(ostream &out=cout)
static std::string className
Name of the class.
Definition: DriverArpack.h:61
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:88
int m_nequ
number of equations
Definition: Driver.h:97
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215
int m_nits
Number of vectors to test.
Definition: DriverArnoldi.h:56