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 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Arnoldi solver using Arpack
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
36 
37 using namespace std;
38 
39 namespace Nektar
40 {
41 namespace SolverUtils
42 {
43 
44 std::string DriverArpack::arpackProblemTypeLookupIds[6] = {
45  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType",
46  "LargestReal", 0),
47  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType",
48  "SmallestReal", 1),
49  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType",
50  "LargestImag", 2),
51  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType",
52  "SmallestImag", 3),
53  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType",
54  "LargestMag", 4),
55  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType",
56  "SmallestMag", 5),
57 };
58 std::string DriverArpack::arpackProblemTypeDefault =
59  LibUtilities::SessionReader::RegisterDefaultSolverInfo("ArpackProblemType",
60  "LargestMag");
61 std::string DriverArpack::driverLookupId =
62  LibUtilities::SessionReader::RegisterEnumValue("Driver", "Arpack", 0);
63 
64 std::string DriverArpack::className =
65  GetDriverFactory().RegisterCreatorFunction("Arpack", DriverArpack::create);
66 
67 std::string DriverArpack::ArpackProblemTypeTrans[6] = {"LR", "SR", "LI",
68  "SI", "LM", "SM"};
69 
70 /**
71  *
72  */
73 DriverArpack::DriverArpack(const LibUtilities::SessionReaderSharedPtr pSession,
75  : DriverArnoldi(pSession, pGraph)
76 {
77 }
78 
79 /**
80  *
81  */
83 {
84 }
85 
86 /**
87  *
88  */
89 void DriverArpack::v_InitObject(ostream &out)
90 {
92 
93  // Initialisation of Arnoldi parameters
94  m_maxn = 1000000; // Maximum size of the problem
95  m_maxnev = 200; // maximum number of eigenvalues requested
96  m_maxncv = 500; // Largest number of basis vector used in Implicitly
97  // Restarted Arnoldi
98 
99  // Error alerts
100  ASSERTL0(m_nvec <= m_maxnev, "NEV is greater than MAXNEV");
101  ASSERTL0(m_kdim <= m_maxncv, "NEV is greater than MAXNEV");
102  ASSERTL0(2 <= m_kdim - m_nvec, "NCV-NEV is less than 2");
103 
104  ASSERTL0(m_comm->GetSize() == 1,
105  "..ARPACK is not currently set-up for parallel execution...\n");
106  m_equ[0]->PrintSummary(out);
107 
108  ASSERTL0(m_comm->GetSize() == 1,
109  "ARPACK Arnoldi solver does not support execution in parallel.");
110 
111  // Print session parameters
112  out << "\tArnoldi solver type : Arpack" << endl;
113 
114  out << "\tArpack problem type : ";
115  out << ArpackProblemTypeTrans[m_session->GetSolverInfoAsEnum<int>(
116  "ArpackProblemType")]
117  << endl;
119 
120  for (int i = 0; i < m_nequ; ++i)
121  {
122  m_equ[i]->DoInitialise();
123  }
124 
125  // FwdTrans Initial conditions to be in Coefficient Space
126  m_equ[m_nequ - 1]->TransPhysToCoeff();
127 }
128 
129 void DriverArpack::v_Execute(ostream &out)
130 
131 {
132  Array<OneD, NekDouble> tmpworkd;
133 
134  int nq = m_equ[0]
135  ->UpdateFields()[0]
136  ->GetNcoeffs(); // Number of points in the mesh
137  int n = m_nfields * nq; // Number of points in eigenvalue calculation
138  int lworkl = 3 * m_kdim * (m_kdim + 2); // Size of work array
139  int ido; // REVERSE COMMUNICATION parameter. At the first call must be
140  // initialised at 0
141  int info; // do not set initial vector (info=0 random initial vector, info=1
142  // read initial vector from session file)
143 
144  int iparam[11];
145  int ipntr[14];
146 
147  Array<OneD, int> ritzSelect;
150  Array<OneD, NekDouble> workev;
152  NekDouble sigmar, sigmai;
153 
154  Array<OneD, NekDouble> resid(n);
156  Array<OneD, NekDouble> workl(lworkl, 0.0);
157  Array<OneD, NekDouble> workd(3 * n, 0.0);
158 
159  ASSERTL0(n <= m_maxn, "N is greater than MAXN");
160 
161  if (m_session->DefinesFunction("InitialConditions"))
162  {
163  out << "\tInital vector : input file " << endl;
164  info = 1;
166  }
167  else
168  {
169  out << "\tInital vector : random " << endl;
170  info = 0;
171  }
172 
173  char B;
174 
175  iparam[0] = 1; // strategy for shift-invert
176  iparam[1] = 0; // (deprecated)
177  iparam[2] = m_nits; // maximum number of iterations allowed/taken
178  iparam[3] = 1; // blocksize to be used for recurrence
179  iparam[4] = 0; // number of converged ritz eigenvalues
180  iparam[5] = 0; // (deprecated)
181 
182  // Use generalized B matrix for coupled solver.
184  {
185  iparam[6] = 1; // computation mode 1=> matrix-vector prod
186  B = 'I';
187  }
188  else
189  {
190  iparam[6] = 3; // computation mode 1=> matrix-vector prod
191  B = 'G';
192  }
193 #if 0
194  if((fabs(m_realShift) > NekConstants::kNekZeroTol)|| // use shift if m_realShift > 1e-12
196  {
197  iparam[6] = 3; // This was 3 need to know what to set it to
198  B = 'G';
199  }
200  else
201  {
202  iparam[6] = 1; // computation mode 1=> matrix-vector prod
203  B = 'I';
204  }
205 #endif
206  iparam[7] = 0; // (for shift-invert)
207  iparam[8] = 0; // number of MV operations
208  iparam[9] = 0; // number of BV operations
209  iparam[10] = 0; // number of reorthogonalisation steps
210 
211  int cycle = 0;
212  const char *problem =
213  ArpackProblemTypeTrans[m_session->GetSolverInfoAsEnum<int>(
214  "ArpackProblemType")]
215  .c_str();
216 
217  std::string name = m_session->GetSessionName() + ".evl";
218  ofstream pFile(name.c_str());
219 
220  ido = 0; // At the first call must be initialisedat 0
221 
222  while (ido != 99) // ido==-1 || ido==1 || ido==0)
223  {
224  // Routine for eigenvalue evaluation for non-symmetric operators
225  Arpack::Dnaupd(ido, &B, // B='I' for std eval problem
226  n, problem, m_nvec, m_evtol, &resid[0], m_kdim, &v[0], n,
227  iparam, ipntr, &workd[0], &workl[0], lworkl, info);
228 
229  // Plotting of real and imaginary part of the
230  // eigenvalues from workl
231  out << "\rIteration " << cycle << ", output: " << info
232  << ", ido=" << ido << " " << std::flush;
233 
234  if (!((cycle - 1) % m_kdim) && (cycle > m_kdim) && (ido != 2))
235  {
236  pFile << "Krylov spectrum at iteration: " << cycle << endl;
237 
239  {
240  pFile << "EV Magnitude Angle Growth Frequency "
241  "Residual"
242  << endl;
243  }
244  else
245  {
246  pFile << "EV Real Imaginary inverse real inverse "
247  "imag Residual"
248  << endl;
249  }
250 
251  out << endl;
252  for (int k = 0; k < m_kdim; ++k)
253  {
254  // write m_kdim eigs to screen
255  WriteEvs(pFile, k, workl[ipntr[5] - 1 + k],
256  workl[ipntr[6] - 1 + k]);
257  }
258  }
259 
260  if (ido == 99)
261  break;
262 
263  switch (ido)
264  {
265  case -1:
266  case 1: // Note that ido=1 we are using input x
267  // (workd[inptr[0]-1]) rather than Mx as
268  // recommended in manual since it is not
269  // possible to impose forcing directly.
270  CopyArnoldiArrayToField(tmpworkd = workd + (ipntr[0] - 1));
271 
272  m_equ[0]->TransCoeffToPhys();
273 
274  m_equ[0]->DoSolve();
276  {
277  // start Adjoint with latest fields of direct
278  CopyFwdToAdj();
279 
280  m_equ[1]->TransCoeffToPhys();
281  m_equ[1]->DoSolve();
282  }
283 
284  if (!(cycle % m_infosteps))
285  {
286  out << endl;
287  m_equ[0]->Output();
288  }
289 
290  // operated fields are copied into workd[inptr[1]-1]
291  CopyFieldToArnoldiArray(tmpworkd = workd + (ipntr[1] - 1));
292 
293  cycle++;
294  break;
295  case 2: // provide y = M x (bwd trans and iproduct);
296  {
297  // workd[inptr[0]-1] copied into operator fields
298  CopyArnoldiArrayToField(tmpworkd = workd + (ipntr[0] - 1));
299 
300  m_equ[0]->TransCoeffToPhys();
301 
303  m_equ[0]->UpdateFields();
304  for (int i = 0; i < fields.size(); ++i)
305  {
306  fields[i]->IProductWRTBase(fields[i]->GetPhys(),
307  fields[i]->UpdateCoeffs());
308  }
309 
310  // operated fields are copied into workd[inptr[1]-1]
311  CopyFieldToArnoldiArray(tmpworkd = workd + (ipntr[1] - 1));
312  break;
313  }
314  default:
315  ASSERTL0(false, "Unexpected reverse communication request.");
316  }
317  }
318 
319  out << endl << "Converged in " << iparam[8] << " iterations" << endl;
320 
321  ASSERTL0(info >= 0, " Error with Dnaupd");
322 
323  ritzSelect = Array<OneD, int>(m_kdim, 0);
324  dr = Array<OneD, NekDouble>(m_nvec + 1, 0.0);
325  di = Array<OneD, NekDouble>(m_nvec + 1, 0.0);
326  workev = Array<OneD, NekDouble>(3 * m_kdim);
327  z = Array<OneD, NekDouble>(n * (m_nvec + 1));
328 
329  if (m_negatedOp)
330  {
331  sigmar = -m_realShift;
332  }
333  else
334  {
335  sigmar = m_realShift;
336  }
337 
338  // Do not pass imaginary shift to Arpack since we have not
339  // used a Fortran complex number format and so processing
340  // is mucked up. Need to do some processing afterwards.
341  sigmai = 0;
342 
343  // Setting 'A', Ritz vectors are computed. 'S' for Shur vectors
344  Arpack::Dneupd(1, "A", ritzSelect.get(), dr.get(), di.get(), z.get(), n,
345  sigmar, sigmai, workev.get(), &B, n, problem, m_nvec,
346  m_evtol, resid.get(), m_kdim, v.get(), n, iparam, ipntr,
347  workd.get(), workl.get(), lworkl, info);
348 
349  ASSERTL0(info == 0, " Error with Dneupd");
350 
351  int nconv = iparam[4];
352 
353  // Subtract off complex shift if it exists
354  if (m_negatedOp)
355  {
356  Vmath::Sadd(nconv, m_imagShift, di, 1, di, 1);
357  }
358  else
359  {
360  Vmath::Sadd(nconv, -m_imagShift, di, 1, di, 1);
361  }
362 
363  WARNINGL0(m_imagShift == 0, "Complex Shift applied. "
364  "Need to implement Ritz re-evaluation of"
365  "eigenvalue. Only one half of complex "
366  "value will be correct");
367 
369  m_equ[0]->UpdateFields();
370 
371  out << "Converged Eigenvalues: " << nconv << endl;
372  pFile << "Converged Eigenvalues: " << nconv << endl;
373 
375  {
376  out << " Magnitude Angle Growth Frequency" << endl;
377  pFile << " Magnitude Angle Growth Frequency"
378  << endl;
379  for (int i = 0; i < nconv; ++i)
380  {
381  WriteEvs(out, i, dr[i], di[i]);
382  WriteEvs(pFile, i, dr[i], di[i]);
383 
384  std::string file = m_session->GetSessionName() + "_eig_" +
385  boost::lexical_cast<std::string>(i) + ".fld";
386  WriteFld(file, z + i * n);
387  }
388  }
389  else
390  {
391  out << " Real Imaginary " << endl;
392  pFile << " Real Imaginary " << endl;
393  for (int i = 0; i < nconv; ++i)
394  {
395  WriteEvs(out, i, dr[i], di[i], NekConstants::kNekUnsetDouble,
396  false);
397  WriteEvs(pFile, i, dr[i], di[i], NekConstants::kNekUnsetDouble,
398  false);
399 
400  std::string file = m_session->GetSessionName() + "_eig_" +
401  boost::lexical_cast<std::string>(i) + ".fld";
402  WriteFld(file, z + i * n);
403  }
404  }
405 
406  m_real_evl = dr;
407  m_imag_evl = di;
408 
409  pFile.close();
410 
411  for (int j = 0; j < m_equ[0]->GetNvariables(); ++j)
412  {
413  NekDouble vL2Error = m_equ[0]->L2Error(j, false);
414  NekDouble vLinfError = m_equ[0]->LinfError(j);
415  if (m_comm->GetRank() == 0)
416  {
417  out << "L 2 error (variable " << m_equ[0]->GetVariable(j)
418  << ") : " << vL2Error << endl;
419  out << "L inf error (variable " << m_equ[0]->GetVariable(j)
420  << ") : " << vLinfError << endl;
421  }
422  }
423 }
424 
425 } // namespace SolverUtils
426 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:222
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
Base class for the development of solvers.
Definition: DriverArnoldi.h:47
void CopyFwdToAdj()
Copy the forward field to the adjoint system in transient growth calculations.
virtual void v_InitObject(std::ostream &out=std::cout) override
int m_infosteps
underlying operator is time stepping
Definition: DriverArnoldi.h:67
void WriteFld(std::string file, std::vector< Array< OneD, NekDouble >> coeffs)
Write coefficients to file.
void CopyFieldToArnoldiArray(Array< OneD, NekDouble > &array)
Copy fields to Arnoldi storage.
int m_nvec
Dimension of Krylov subspace.
Definition: DriverArnoldi.h:61
bool m_timeSteppingAlgorithm
Period of time stepping algorithm.
Definition: DriverArnoldi.h:65
int m_nits
Number of vectors to test.
Definition: DriverArnoldi.h:62
Array< OneD, NekDouble > m_imag_evl
Definition: DriverArnoldi.h:75
void CopyArnoldiArrayToField(Array< OneD, NekDouble > &array)
Copy Arnoldi storage to fields.
NekDouble m_evtol
Maxmum number of iterations.
Definition: DriverArnoldi.h:63
int m_nfields
interval to dump information if required.
Definition: DriverArnoldi.h:69
SOLVER_UTILS_EXPORT void ArnoldiSummary(std::ostream &out)
Array< OneD, NekDouble > m_real_evl
Operator in solve call is negated.
Definition: DriverArnoldi.h:74
void WriteEvs(std::ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble, bool DumpInverse=true)
virtual void v_InitObject(std::ostream &out=std::cout) override
Virtual function for initialisation implementation.
virtual ~DriverArpack()
Destructor.
virtual void v_Execute(std::ostream &out=std::cout) override
Virtual function for solve implementation.
static std::string ArpackProblemTypeTrans[]
Definition: DriverArpack.h:91
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:88
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:85
enum EvolutionOperatorType m_EvolutionOperator
Evolution Operator.
Definition: Driver.h:103
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:97
int m_nequ
number of equations
Definition: Driver.h:100
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:116
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:131
std::shared_ptr< SessionReader > SessionReaderSharedPtr
static const NekDouble kNekUnsetDouble
static const NekDouble kNekZeroTol
DriverFactory & GetDriverFactory()
Definition: Driver.cpp:64
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add scalar y = alpha + x.
Definition: Vmath.cpp:384