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 {
42 
44  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","LargestReal" ,0),
45  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","SmallestReal" ,1),
46  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","LargestImag" ,2),
47  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","SmallestImag" ,3),
48  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","LargestMag" ,4),
49  LibUtilities::SessionReader::RegisterEnumValue("ArpackProblemType","SmallestMag" ,5),
50 };
53 
55 
57 { "LR", "SR", "LI", "SI", "LM", "SM" };
58 
59 
60 /**
61  *
62  */
64  : DriverArnoldi(pSession)
65 {
66 }
67 
68 
69 /**
70  *
71  */
73 {
74 }
75 
76 /**
77  *
78  */
79 void DriverArpack::v_InitObject(ostream &out)
80 {
82 
83  //Initialisation of Arnoldi parameters
84  m_maxn = 1000000; // Maximum size of the problem
85  m_maxnev = 200; // maximum number of eigenvalues requested
86  m_maxncv = 500; // Largest number of basis vector used in Implicitly Restarted Arnoldi
87 
88  // Error alerts
89  ASSERTL0(m_nvec <= m_maxnev,"NEV is greater than MAXNEV");
90  ASSERTL0(m_kdim <= m_maxncv,"NEV is greater than MAXNEV");
91  ASSERTL0(2 <= m_kdim-m_nvec,"NCV-NEV is less than 2");
92 
93  m_equ[0]->PrintSummary(out);
94 
95  // Print session parameters
96  out << "\tArnoldi solver type : Arpack" << endl;
97 
98  out << "\tArpack problem type : ";
99  out << ArpackProblemTypeTrans[m_session->GetSolverInfoAsEnum<int>("ArpackProblemType")] << endl;
101 
102  m_equ[m_nequ - 1]->DoInitialise();
103 
104  // FwdTrans Initial conditions to be in Coefficient Space
105  m_equ[m_nequ-1] ->TransPhysToCoeff();
106 
107 }
108 
109 
110 void DriverArpack::v_Execute(ostream &out)
111 
112 {
113  Array<OneD, NekDouble> tmpworkd;
114 
115  int nq = m_equ[0]->UpdateFields()[0]->GetNcoeffs(); // Number of points in the mesh
116  int n = m_nfields*nq; // Number of points in eigenvalue calculation
117  int lworkl = 3*m_kdim*(m_kdim+2); // Size of work array
118  int ido ; //REVERSE COMMUNICATION parameter. At the first call must be initialised at 0
119  int info; // do not set initial vector (info=0 random initial vector, info=1 read initial vector from session file)
120 
121  int iparam[11];
122  int ipntr[14];
123 
124  Array<OneD, int> ritzSelect;
127  Array<OneD, NekDouble> workev;
129  NekDouble sigmar, sigmai;
130 
131  Array<OneD, NekDouble> resid(n);
133  Array<OneD, NekDouble> workl(lworkl, 0.0);
134  Array<OneD, NekDouble> workd(3*n, 0.0);
135 
136  ASSERTL0(n <= m_maxn, "N is greater than MAXN");
137 
138  if(m_session->DefinesFunction("InitialConditions"))
139  {
140  out << "\tInital vector : input file " << endl;
141  info = 1;
143  }
144  else
145  {
146  out << "\tInital vector : random " << endl;
147  info = 0;
148  }
149 
150  char B;
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 
159  // Use generalized B matrix for coupled solver.
161  {
162  iparam[6] = 1; // computation mode 1=> matrix-vector prod
163  B = 'I';
164  }
165  else {
166  iparam[6] = 3; // computation mode 1=> matrix-vector prod
167  B = 'G';
168  }
169 #if 0
170  if((fabs(m_realShift) > NekConstants::kNekZeroTol)|| // use shift if m_realShift > 1e-12
172  {
173  iparam[6] = 3; // This was 3 need to know what to set it to
174  B = 'G';
175  }
176  else
177  {
178  iparam[6] = 1; // computation mode 1=> matrix-vector prod
179  B = 'I';
180  }
181 #endif
182  iparam[7] = 0; // (for shift-invert)
183  iparam[8] = 0; // number of MV operations
184  iparam[9] = 0; // number of BV operations
185  iparam[10]= 0; // number of reorthogonalisation steps
186 
187  int cycle = 0;
188  const char* problem = ArpackProblemTypeTrans[m_session->GetSolverInfoAsEnum<int>("ArpackProblemType")].c_str();
189 
190  std::string name = m_session->GetSessionName() + ".evl";
191  ofstream pFile(name.c_str());
192 
193  ido = 0; //At the first call must be initialisedat 0
194 
195  while(ido != 99)//ido==-1 || ido==1 || ido==0)
196  {
197  //Routine for eigenvalue evaluation for non-symmetric operators
198  Arpack::Dnaupd( ido, &B, // B='I' for std eval problem
199  n, problem, m_nvec,
200  m_evtol, &resid[0], m_kdim,
201  &v[0], n, iparam, ipntr, &workd[0],
202  &workl[0], lworkl, info);
203 
204  //Plotting of real and imaginary part of the
205  //eigenvalues from workl
206  out << "\rIteration " << cycle << ", output: " << info << ", ido=" << ido << " " << std::flush;
207 
208  if(!((cycle-1)%m_kdim)&&(cycle> m_kdim)&&(ido!=2))
209  {
210  pFile << "Krylov spectrum at iteration: " << cycle << endl;
211 
213  {
214  pFile << "EV Magnitude Angle Growth Frequency Residual" << endl;
215  }
216  else
217  {
218  pFile << "EV Real Imaginary inverse real inverse imag Residual" << endl;
219  }
220 
221  out << endl;
222  for(int k = 0; k < m_kdim; ++k)
223  {
224  // write m_kdim eigs to screen
225  WriteEvs(pFile,k, workl[ipntr[5]-1+k],workl[ipntr[6]-1+k]);
226  }
227  }
228 
229  if (ido == 99) break;
230 
231  switch(ido)
232  {
233  case -1:
234  case 1: // Note that ido=1 we are using input x
235  // (workd[inptr[0]-1]) rather than Mx as
236  // recommended in manual since it is not
237  // possible to impose forcing directly.
238  CopyArnoldiArrayToField(tmpworkd = workd + (ipntr[0]-1));
239 
240  m_equ[0]->TransCoeffToPhys();
241 
242  m_equ[0]->DoSolve();
244  {
245  //start Adjoint with latest fields of direct
246  CopyFwdToAdj();
247 
248  m_equ[1]->TransCoeffToPhys();
249  m_equ[1]->DoSolve();
250  }
251 
252  if(!(cycle%m_infosteps))
253  {
254  out << endl;
255  m_equ[0]->Output();
256  }
257 
258  // operated fields are copied into workd[inptr[1]-1]
259  CopyFieldToArnoldiArray(tmpworkd = workd + (ipntr[1]-1));
260 
261  cycle++;
262  break;
263  case 2: // provide y = M x (bwd trans and iproduct);
264  {
265  //workd[inptr[0]-1] copied into operator fields
266  CopyArnoldiArrayToField(tmpworkd = workd + (ipntr[0]-1));
267 
268  m_equ[0]->TransCoeffToPhys();
269 
270  Array<OneD, MultiRegions::ExpListSharedPtr> fields = m_equ[0]->UpdateFields();
271  for (int i = 0; i < fields.num_elements(); ++i)
272  {
273  fields[i]->IProductWRTBase(fields[i]->GetPhys(),
274  fields[i]->UpdateCoeffs());
275  }
276 
277  // operated fields are copied into workd[inptr[1]-1]
278  CopyFieldToArnoldiArray(tmpworkd = workd + (ipntr[1]-1));
279  break;
280  }
281  default:
282  ASSERTL0(false, "Unexpected reverse communication request.");
283  }
284 
285  }
286 
287  out<< endl << "Converged in " << iparam[8] << " iterations" << endl;
288 
289  ASSERTL0(info >= 0," Error with Dnaupd");
290 
291  ritzSelect = Array<OneD, int> (m_kdim,0);
292  dr = Array<OneD, NekDouble> (m_nvec+1,0.0);
293  di = Array<OneD, NekDouble> (m_nvec+1,0.0);
294  workev = Array<OneD, NekDouble> (3*m_kdim);
295  z = Array<OneD, NekDouble> (n*(m_nvec+1));
296 
297  if(m_negatedOp)
298  {
299  sigmar = -m_realShift;
300  }
301  else
302  {
303  sigmar = m_realShift;
304  }
305 
306  // Do not pass imaginary shift to Arpack since we have not
307  // used a Fortran complex number format and so processing
308  // is mucked up. Need to do some processing afterwards.
309  sigmai = 0;
310 
311  //Setting 'A', Ritz vectors are computed. 'S' for Shur vectors
312  Arpack::Dneupd(1, "A", ritzSelect.get(), dr.get(), di.get(),
313  z.get(), n, sigmar, sigmai, workev.get(), &B,
314  n, problem, m_nvec, m_evtol, resid.get(), m_kdim,
315  v.get(), n, iparam, ipntr, workd.get(),
316  workl.get(),lworkl,info);
317 
318  ASSERTL0(info == 0, " Error with Dneupd");
319 
320  int nconv=iparam[4];
321 
322  // Subtract off complex shift if it exists
323  if(m_negatedOp)
324  {
325  Vmath::Sadd(nconv,m_imagShift,di,1,di,1);
326  }
327  else
328  {
329  Vmath::Sadd(nconv,-m_imagShift,di,1,di,1);
330  }
331 
332  WARNINGL0(m_imagShift == 0,"Complex Shift applied. "
333  "Need to implement Ritz re-evaluation of"
334  "eigenvalue. Only one half of complex "
335  "value will be correct");
336 
337 
338  Array<OneD, MultiRegions::ExpListSharedPtr> fields = m_equ[0]->UpdateFields();
339 
340  out << "Converged Eigenvalues: " << nconv << endl;
341  pFile << "Converged Eigenvalues: " << nconv << endl;
342 
344  {
345  out << " Magnitude Angle Growth Frequency" << endl;
346  pFile << " Magnitude Angle Growth Frequency" << endl;
347  for(int i= 0; i< nconv; ++i)
348  {
349  WriteEvs(out,i,dr[i],di[i]);
350  WriteEvs(pFile,i,dr[i],di[i]);
351 
352  std::string file = m_session->GetSessionName() + "_eig_"
353  + boost::lexical_cast<std::string>(i)
354  + ".fld";
355  WriteFld(file,z + i*n);
356  }
357  }
358  else
359  {
360  out << " Real Imaginary " << endl;
361  pFile << " Real Imaginary " << endl;
362  for(int i= 0; i< nconv; ++i)
363  {
364  WriteEvs(out,i,dr[i],di[i],
366  WriteEvs(pFile,i,dr[i],di[i],
368 
369  std::string file = m_session->GetSessionName() + "_eig_"
370  + boost::lexical_cast<std::string>(i)
371  + ".fld";
372  WriteFld(file,z + i*n);
373  }
374  }
375 
376  m_real_evl = dr;
377  m_imag_evl = di;
378 
379  pFile.close();
380 
381  for(int j = 0; j < m_equ[0]->GetNvariables(); ++j)
382  {
383  NekDouble vL2Error = m_equ[0]->L2Error(j,false);
384  NekDouble vLinfError = m_equ[0]->LinfError(j);
385  if (m_comm->GetRank() == 0)
386  {
387  out << "L 2 error (variable " << m_equ[0]->GetVariable(j) << ") : " << vL2Error << endl;
388  out << "L inf error (variable " << m_equ[0]->GetVariable(j) << ") : " << vLinfError << endl;
389  }
390  }
391 }
392 
393 }
394 }
Array< OneD, NekDouble > m_imag_evl
Definition: DriverArnoldi.h:70
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:161
void CopyFwdToAdj()
Copy the forward field to the adjoint system in transient growth calculations.
NekDouble m_evtol
Maxmum number of iterations.
Definition: DriverArnoldi.h:58
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)
Write coefficients to file.
static std::string arpackProblemTypeLookupIds[]
Definition: DriverArpack.h:83
DriverArpack(const LibUtilities::SessionReaderSharedPtr pSession)
Constructor.
static const NekDouble kNekZeroTol
int m_nvec
Dimension of Krylov subspace.
Definition: DriverArnoldi.h:56
Array< OneD, NekDouble > m_real_evl
Operator in solve call is negated.
Definition: DriverArnoldi.h:69
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:167
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
static const NekDouble kNekUnsetDouble
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.cpp:301
int m_nfields
interval to dump information if required.
Definition: DriverArnoldi.h:64
bool m_timeSteppingAlgorithm
Period of time stepping algorithm.
Definition: DriverArnoldi.h:60
void WriteEvs(ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble, bool DumpInverse=true)
static std::string ArpackProblemTypeTrans[]
Definition: DriverArpack.h:88
static std::string driverLookupId
Definition: DriverArpack.h:85
static std::string arpackProblemTypeDefault
Definition: DriverArpack.h:84
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:94
int m_infosteps
underlying operator is time stepping
Definition: DriverArnoldi.h:62
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:55
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:47
virtual void v_InitObject(ostream &out=cout)
static std::string className
Name of the class.
Definition: DriverArpack.h:62
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:57