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