Nektar++
DriverModifiedArnoldi.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File DriverModifiedArnoldi.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: Driver for eigenvalue analysis using the modified Arnoldi
32 // method.
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include <iomanip>
37 
38 #include <boost/core/ignore_unused.hpp>
39 
41 
42 using namespace std;
43 
44 namespace Nektar
45 {
46 namespace SolverUtils
47 {
48 
49 string DriverModifiedArnoldi::className =
50  GetDriverFactory().RegisterCreatorFunction("ModifiedArnoldi",
51  DriverModifiedArnoldi::create);
52 string DriverModifiedArnoldi::driverLookupId =
53  LibUtilities::SessionReader::RegisterEnumValue("Driver",
54  "ModifiedArnoldi",0);
55 
56 /**
57  *
58  */
59 DriverModifiedArnoldi::DriverModifiedArnoldi(
62  : DriverArnoldi(pSession, pGraph)
63 {
64 }
65 
66 
67 /**
68  *
69  */
71 {
72 }
73 
74 
75 /**
76  *
77  */
79 {
81 
82  m_equ[0]->PrintSummary(out);
83 
84  // Print session parameters
85  if (m_comm->GetRank() == 0)
86  {
87  out << "\tArnoldi solver type : Modified Arnoldi" << endl;
88  }
89 
91 
92  for( int i = 0; i < m_nequ; ++i)
93  {
94  m_equ[i]->DoInitialise();
95  }
96 
97  //FwdTrans Initial conditions to be in Coefficient Space
98  m_equ[m_nequ-1] ->TransPhysToCoeff();
99 
100 }
101 
102 
103 /**
104  *
105  */
107 {
108  int i = 0;
109  int j = 0;
110  int nq = m_equ[0]->UpdateFields()[0]->GetNcoeffs();
111  int ntot = m_nfields*nq;
112  int converged = 0;
113  NekDouble resnorm = 0.0;
114  ofstream evlout;
115  std::string evlFile = m_session->GetSessionName() + ".evl";
116 
117  if (m_comm->GetRank() == 0)
118  {
119  evlout.open(evlFile.c_str());
120  }
121 
122  // Allocate memory
127 
134  for (i = 0; i < m_kdim + 1; ++i)
135  {
136  Kseq[i] = Array<OneD, NekDouble>(ntot, 0.0);
137  if (m_useMask)
138  {
139  Kseqcopy[i] = Array<OneD, NekDouble>(ntot, 0.0);
140  }
141  else
142  {
143  Kseqcopy[i] = Kseq[i];
144  }
145  Tseq[i] = Array<OneD, NekDouble>(ntot, 0.0);
146  }
147 
148  // Copy starting vector into second sequence element (temporary).
149  if(m_session->DefinesFunction("InitialConditions"))
150  {
151  if (m_comm->GetRank() == 0)
152  {
153  out << "\tInital vector : specified in input file " << endl;
154  }
155  m_equ[0]->SetInitialConditions(0.0,false);
156 
157  CopyFieldToArnoldiArray(Kseq[1]);
158  }
159  else
160  {
161  if (m_comm->GetRank() == 0)
162  {
163  out << "\tInital vector : random " << endl;
164  }
165 
166  NekDouble eps=0.0001;
167  Vmath::FillWhiteNoise(ntot, eps , &Kseq[1][0], 1);
168  if (m_useMask)
169  {
170  Vmath::Vmul(ntot, Kseq[1], 1, GetMaskCoeff(), 1, Kseq[1], 1);
171  }
172  }
173 
174  // Perform one iteration to enforce boundary conditions.
175  // Set this as the initial value in the sequence.
176  EV_update(Kseq[1], Kseq[0]);
177  if (m_comm->GetRank() == 0)
178  {
179  out << "Iteration: " << 0 << endl;
180  }
181 
182  // Normalise first vector in sequence
183  if (m_useMask)
184  {
185  Vmath::Vmul(ntot, Kseq[0], 1, GetMaskCoeff(), 1, Kseqcopy[0], 1);
186  }
187  alpha[0] = Blas::Ddot(ntot, &Kseqcopy[0][0], 1, &Kseqcopy[0][0], 1);
188  m_comm->AllReduce(alpha[0], Nektar::LibUtilities::ReduceSum);
189  alpha[0] = std::sqrt(alpha[0]);
190  Vmath::Smul(ntot, 1.0/alpha[0], Kseq[0], 1, Kseq[0], 1);
191 
192  // Fill initial krylov sequence
193  NekDouble resid0;
194  for (i = 1; !converged && i <= m_kdim; ++i)
195  {
196  // Compute next vector
197  EV_update(Kseq[i-1], Kseq[i]);
198 
199  // Normalise
200  if (m_useMask)
201  {
202  Vmath::Vmul(ntot, Kseq[i], 1, GetMaskCoeff(), 1, Kseqcopy[i], 1);
203  }
204  alpha[i] = Blas::Ddot(ntot, &Kseqcopy[i][0], 1, &Kseqcopy[i][0], 1);
205  m_comm->AllReduce(alpha[i], Nektar::LibUtilities::ReduceSum);
206  alpha[i] = std::sqrt(alpha[i]);
207 
208  //alpha[i] = std::sqrt(alpha[i]);
209  Vmath::Smul(ntot, 1.0/alpha[i], Kseq[i], 1, Kseq[i], 1);
210 
211  // Copy Krylov sequence into temporary storage
212  for (int k = 0; k < i + 1; ++k)
213  {
214  if (m_useMask)
215  {
216  Vmath::Vmul(ntot, Kseq[k], 1, GetMaskCoeff(), 1, Tseq[k], 1);
217  Vmath::Vcopy(ntot, Kseq[k], 1, Kseqcopy[k], 1);
218  }
219  else
220  {
221  Vmath::Vcopy(ntot, Kseq[k], 1, Tseq[k], 1);;
222  }
223  }
224 
225  // Generate Hessenberg matrix and compute eigenvalues of it.
226  EV_small(Tseq, Kseqcopy, ntot, alpha, i, zvec, wr, wi, resnorm);
227 
228  // Test for convergence.
229  converged = EV_test(i, i, zvec, wr, wi, resnorm,
230  std::min(i, m_nvec), evlout, resid0);
231  if ( i >= m_nvec)
232  {
233  converged = max (converged, 0);
234  }
235  else
236  {
237  converged = 0;
238  }
239 
240  if (m_comm->GetRank() == 0)
241  {
242  out << "Iteration: " << i << " (residual : " << resid0
243  << ")" <<endl;
244  }
245  }
246 
247  // Continue with full sequence
248  if (!converged)
249  {
250  for (i = m_kdim + 1; !converged && i <= m_nits; ++i)
251  {
252  // Shift all the vectors in the sequence.
253  // First vector is removed.
254  for (int j = 1; j <= m_kdim; ++j)
255  {
256  alpha[j-1] = alpha[j];
257  Vmath::Vcopy(ntot, Kseq[j], 1, Kseq[j-1], 1);
258  }
259 
260  // Compute next vector
261  EV_update(Kseq[m_kdim - 1], Kseq[m_kdim]);
262 
263  // Compute new scale factor
264  if (m_useMask)
265  {
266  Vmath::Vmul(ntot, Kseq[m_kdim], 1, GetMaskCoeff(), 1,
267  Kseqcopy[m_kdim], 1);
268  }
269  alpha[m_kdim] = Blas::Ddot(ntot, &Kseqcopy[m_kdim][0], 1,
270  &Kseqcopy[m_kdim][0], 1);
271  m_comm->AllReduce(alpha[m_kdim], Nektar::LibUtilities::ReduceSum);
272  alpha[m_kdim] = std::sqrt(alpha[m_kdim]);
273  Vmath::Smul(ntot, 1.0/alpha[m_kdim], Kseq[m_kdim], 1,
274  Kseq[m_kdim], 1);
275 
276  // Copy Krylov sequence into temporary storage
277  for (int k = 0; k < m_kdim + 1; ++k)
278  {
279  if (m_useMask)
280  {
281  Vmath::Vmul(ntot, Kseq[k], 1, GetMaskCoeff(), 1, Tseq[k], 1);
282  Vmath::Vcopy(ntot, Kseq[k], 1, Kseqcopy[k], 1);
283  }
284  else
285  {
286  Vmath::Vcopy(ntot, Kseq[k], 1, Tseq[k], 1);;
287  }
288  }
289 
290  // Generate Hessenberg matrix and compute eigenvalues of it
291  EV_small(Tseq, Kseqcopy, ntot, alpha, m_kdim, zvec, wr, wi, resnorm);
292 
293  // Test for convergence.
294  converged = EV_test(i, m_kdim, zvec, wr, wi, resnorm,
295  m_nvec, evlout, resid0);
296 
297  if (m_comm->GetRank() == 0)
298  {
299  out << "Iteration: " << i << " (residual : "
300  << resid0 << ")" <<endl;
301  }
302  }
303  }
304 
305  m_equ[0]->Output();
306 
307  // Evaluate and output computation time and solution accuracy.
308  // The specific format of the error output is essential for the
309  // regression tests to work.
310  // Evaluate L2 Error
311  for(j = 0; j < m_equ[0]->GetNvariables(); ++j)
312  {
313  NekDouble vL2Error = m_equ[0]->L2Error(j,false);
314  NekDouble vLinfError = m_equ[0]->LinfError(j);
315  if (m_comm->GetRank() == 0)
316  {
317  out << "L 2 error (variable " << m_equ[0]->GetVariable(j)
318  << ") : " << vL2Error << endl;
319  out << "L inf error (variable " << m_equ[0]->GetVariable(j)
320  << ") : " << vLinfError << endl;
321  }
322  }
323 
324  // Process eigenvectors and write out.
325  m_nvec = converged;
326  EV_post(Tseq, Kseqcopy, ntot, min(--i, m_kdim), m_nvec, zvec, wr, wi,
327  converged);
328 
329  WARNINGL0(m_imagShift == 0,"Complex Shift applied. "
330  "Need to implement Ritz re-evaluation of"
331  "eigenvalue. Only one half of complex "
332  "value will be correct");
333 
334  // store eigenvalues so they can be accessed from driver class
335  m_real_evl = wr;
336  m_imag_evl = wi;
337 
338  // Close the runtime info file.
339  if (m_comm->GetRank() == 0)
340  {
341  evlout.close();
342  }
343 }
344 
345 
346 /**
347  *
348  */
352 {
353  // Copy starting vector into first sequence element.
355  m_equ[0]->TransCoeffToPhys();
356 
357  m_equ[0]->SetTime(0.);
358  m_equ[0]->DoSolve();
359 
361  {
363  fields = m_equ[0]->UpdateFields();
364 
365  //start Adjoint with latest fields of direct
366  CopyFwdToAdj();
367  m_equ[1]->TransCoeffToPhys();
368 
369  m_equ[1]->SetTime(0.);
370  m_equ[1]->DoSolve();
371  }
372 
373  // Copy starting vector into first sequence element.
375 }
376 
377 
378 /**
379  * Computes the Ritz eigenvalues and eigenvectors of the Hessenberg matrix
380  * constructed from the Krylov sequence.
381  */
384  Array<OneD, Array<OneD, NekDouble> > &Kseqcopy,
385  const int ntot,
386  const Array<OneD, NekDouble> &alpha,
387  const int kdim,
391  NekDouble &resnorm)
392 {
393  int kdimp = kdim + 1;
394  int lwork = 10*kdim;
395  int ier;
396  Array<OneD, NekDouble> R(kdimp * kdimp, 0.0);
397  Array<OneD, NekDouble> H(kdimp * kdim, 0.0);
398  Array<OneD, NekDouble> rwork(lwork, 0.0);
399 
400  // Modified G-S orthonormalisation
401  for (int i = 0; i < kdimp; ++i)
402  {
403  NekDouble gsc = Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[i][0], 1);
404  m_comm->AllReduce(gsc, Nektar::LibUtilities::ReduceSum);
405  gsc = std::sqrt(gsc);
406  ASSERTL0(gsc != 0.0, "Vectors are linearly independent.");
407 
408  R[i*kdimp+i] = gsc;
409  Vmath::Smul(ntot, 1.0/gsc, Kseq[i], 1, Kseq[i], 1);
410  if (m_useMask)
411  {
412  Vmath::Smul(ntot, 1.0/gsc, Kseqcopy[i], 1, Kseqcopy[i], 1);
413  }
414 
415  for (int j = i + 1; j < kdimp; ++j)
416  {
417  gsc = Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[j][0], 1);
418  m_comm->AllReduce(gsc, Nektar::LibUtilities::ReduceSum);
419  Vmath::Svtvp(ntot, -gsc, Kseq[i], 1, Kseq[j], 1, Kseq[j], 1);
420  if (m_useMask)
421  {
422  Vmath::Svtvp(ntot, -gsc, Kseqcopy[i], 1, Kseqcopy[j], 1, Kseqcopy[j], 1);
423  }
424  R[j*kdimp+i] = gsc;
425  }
426  }
427 
428  // Compute H matrix
429  for (int i = 0; i < kdim; ++i)
430  {
431  for (int j = 0; j < kdim; ++j)
432  {
433  H[j*kdim+i] = alpha[j+1] * R[(j+1)*kdimp+i]
434  - Vmath::Dot(j, &H[0] + i, kdim, &R[0] + j*kdimp, 1);
435  H[j*kdim+i] /= R[j*kdimp+j];
436  }
437  }
438 
439  H[(kdim-1)*kdim+kdim] = alpha[kdim]
440  * std::fabs(R[kdim*kdimp+kdim] / R[(kdim-1)*kdimp + kdim-1]);
441 
442  Lapack::dgeev_('N', 'V', kdim, &H[0], kdim, &wr[0], &wi[0], 0, 1,
443  &zvec[0], kdim, &rwork[0], lwork, ier);
444 
445  ASSERTL0(!ier, "Error with dgeev");
446 
447  resnorm = H[(kdim-1)*kdim + kdim];
448 }
449 
450 
451 /**
452  * Check for convergence of the residuals of the eigenvalues of H.
453  */
455  const int itrn,
456  const int kdim,
460  const NekDouble resnorm,
461  int nvec,
462  ofstream &evlout,
463  NekDouble &resid0)
464 {
465  int idone = 0;
466  // NekDouble period = 0.1;
467 
468  Array<OneD, NekDouble> resid(kdim);
469  for (int i = 0; i < kdim; ++i)
470  {
471  NekDouble tmp = std::sqrt(Vmath::Dot(kdim, &zvec[0] + i*kdim, 1,
472  &zvec[0] + i*kdim, 1));
473  resid[i] = resnorm * std::fabs(zvec[kdim - 1 + i*kdim]) / tmp;
474  if (wi[i] < 0.0)
475  {
476  resid[i-1] = resid[i] = hypot(resid[i-1], resid[i]);
477  }
478  }
479 
480  EV_sort(zvec, wr, wi, resid, kdim);
481 
482  while (nvec <= kdim && resid[nvec-1] < m_evtol)
483  {
484  idone = nvec;
485  ++nvec;
486  }
487  nvec -= (idone > 0);
488 
489  if (m_comm->GetRank() == 0)
490  {
491  evlout << "-- Iteration = " << itrn << ", H(k+1, k) = "
492  << resnorm << endl;
493  evlout.precision(4);
494  evlout.setf(ios::scientific, ios::floatfield);
496  {
497  evlout << " Magnitude Angle Growth "
498  << "Frequency Residual" << endl;
499  }
500  else
501  {
502  evlout << " Real Imaginary inverse real "
503  << "inverse imag Residual" << endl;
504  }
505 
506  for (int i = 0; i < kdim; i++)
507  {
508  WriteEvs(evlout,i,wr[i],wi[i],resid[i]);
509  }
510  }
511 
512  resid0 = resid[nvec-1];
513  return idone;
514 }
515 
516 
517 /**
518  * Sorts the computed eigenvalues by smallest residual first.
519  */
525  const int dim)
526 {
527  Array<OneD, NekDouble> z_tmp(dim,0.0);
528  NekDouble wr_tmp, wi_tmp, te_tmp;
529  for (int j = 1; j < dim; ++j)
530  {
531  wr_tmp = wr[j];
532  wi_tmp = wi[j];
533  te_tmp = test[j];
534  Vmath::Vcopy(dim, &evec[0] + j*dim, 1, &z_tmp[0], 1);
535  int i = j - 1;
536  while (i >= 0 && test[i] > te_tmp)
537  {
538  wr[i+1] = wr[i];
539  wi[i+1] = wi[i];
540  test[i+1] = test[i];
541  Vmath::Vcopy(dim, &evec[0] + i*dim, 1, &evec[0] + (i+1)*dim, 1);
542  i--;
543  }
544  wr[i+1] = wr_tmp;
545  wi[i+1] = wi_tmp;
546  test[i+1] = te_tmp;
547  Vmath::Vcopy(dim, &z_tmp[0], 1, &evec[0] + (i+1)*dim, 1);
548  }
549 }
550 
551 
552 /**
553  * Post-process the Ritz eigenvalues/eigenvectors of the matrix H, to compute
554  * estimations of the leading eigenvalues and eigenvectors of the original
555  * matrix.
556  */
560  const int ntot,
561  const int kdim,
562  const int nvec,
566  const int icon)
567 {
568  if (icon == 0)
569  {
570  // Not converged, write final Krylov vector
571  ASSERTL0(false, "Convergence was not achieved within the "
572  "prescribed number of iterations.");
573  }
574  else if (icon < 0)
575  {
576  // Minimum residual reached
577  ASSERTL0(false, "Minimum residual reached.");
578  }
579  else if (icon == nvec)
580  {
581  // Converged, write out eigenvectors
582  EV_big(Tseq, Kseq, ntot, kdim, icon, zvec, wr, wi);
584  = m_equ[0]->UpdateFields();
585 
586  for (int j = 0; j < icon; ++j)
587  {
588  std::string file = m_session->GetSessionName() + "_eig_"
589  + boost::lexical_cast<std::string>(j)
590  + ".fld";
591 
592  if (m_comm->GetRank() == 0)
593  {
594  WriteEvs(cout, j, wr[j], wi[j]);
595  }
596  WriteFld(file,Kseq[j]);
597  if (m_useMask)
598  {
599  std::string fileunmask = m_session->GetSessionName() + "_eig_masked_"
600  + boost::lexical_cast<std::string>(j)
601  + ".fld";
602  WriteFld(fileunmask,Tseq[j]);
603  }
604  }
605  }
606  else
607  {
608  // Not recognised
609  ASSERTL0(false, "Unrecognised value.");
610  }
611 }
612 
613 
614 /**
615  * Compute estimates of the eigenvalues/eigenvectors of the original system.
616  */
620  const int ntot,
621  const int kdim,
622  const int nvec,
626 {
627  boost::ignore_unused(wr);
628 
629  NekDouble wgt, norm;
632  for(int i=0; i<nvec; ++i)
633  {
634  if (m_useMask)
635  {
636  btmp[i] = Array<OneD, NekDouble>(ntot, 0.);
637  etmp[i] = Array<OneD, NekDouble>(ntot, 0.);
638  }
639  else
640  {
641  btmp[i] = evecs[i];
642  }
643  }
644 
645  // Generate big eigenvectors
646  for (int j = 0; j < nvec; ++j)
647  {
648  Vmath::Zero(ntot, btmp[j], 1);
649  if (m_useMask)
650  {
651  Vmath::Zero(ntot, etmp[j], 1);
652  }
653  for (int i = 0; i < kdim; ++i)
654  {
655  wgt = zvec[i + j*kdim];
656  Vmath::Svtvp(ntot, wgt, bvecs[i], 1, btmp[j], 1, btmp[j], 1);
657  if (m_useMask)
658  {
659  Vmath::Svtvp(ntot, wgt, evecs[i], 1, etmp[j], 1, etmp[j], 1);
660  }
661  }
662  }
663 
664  // Normalise the big eigenvectors
665  for (int i = 0; i < nvec; ++i)
666  {
667  if (wi[i] == 0.0) // Real mode
668  {
669  if (m_session->GetComm()->GetRank() == 0)
670  {
671  cout << "eigenvalue " << i << ": real mode" << endl;
672  }
673  norm = Blas::Ddot(ntot, &btmp[i][0], 1, &btmp[i][0], 1);
674  m_comm->AllReduce(norm, Nektar::LibUtilities::ReduceSum);
675  norm = std::sqrt(norm);
676  if (m_useMask)
677  {
678  Vmath::Smul(ntot, 1.0/norm, btmp[i], 1, bvecs[i], 1);
679  Vmath::Smul(ntot, 1.0/norm, etmp[i], 1, evecs[i], 1);
680  }
681  else
682  {
683  Vmath::Smul(ntot, 1.0/norm, btmp[i], 1, evecs[i], 1);
684  }
685  }
686  else
687  {
688  if (m_session->GetComm()->GetRank() == 0)
689  {
690  cout << "eigenvalues " << i << ", " << i+1
691  << ": complex modes" << endl;
692  }
693  norm = Blas::Ddot(ntot, &btmp[i][0], 1, &btmp[i][0], 1);
694  norm += Blas::Ddot(ntot, &btmp[i+1][0], 1, &btmp[i+1][0], 1);
695  m_comm->AllReduce(norm, Nektar::LibUtilities::ReduceSum);
696  norm = std::sqrt(norm);
697 
698  if (m_useMask)
699  {
700  Vmath::Smul(ntot, 1.0/norm, btmp[i], 1, bvecs[i], 1);
701  Vmath::Smul(ntot, 1.0/norm, btmp[i+1], 1, bvecs[i+1], 1);
702  Vmath::Smul(ntot, 1.0/norm, etmp[i], 1, evecs[i], 1);
703  Vmath::Smul(ntot, 1.0/norm, etmp[i+1], 1, evecs[i+1], 1);
704  }
705  else
706  {
707  Vmath::Smul(ntot, 1.0/norm, btmp[i], 1, evecs[i], 1);
708  Vmath::Smul(ntot, 1.0/norm, btmp[i+1], 1, evecs[i+1], 1);
709  }
710 
711  i++;
712  }
713  }
714 }
715 
716 }
717 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:223
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
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.
void WriteFld(std::string file, std::vector< Array< OneD, NekDouble > > coeffs)
Write coefficients to file.
SOLVER_UTILS_EXPORT const Array< OneD, const NekDouble > & GetMaskCoeff() const
void CopyFieldToArnoldiArray(Array< OneD, NekDouble > &array)
Copy fields to Arnoldi storage.
virtual void v_InitObject(std::ostream &out=std::cout)
int m_nvec
Dimension of Krylov subspace.
Definition: DriverArnoldi.h:59
bool m_timeSteppingAlgorithm
Period of time stepping algorithm.
Definition: DriverArnoldi.h:63
int m_nits
Number of vectors to test.
Definition: DriverArnoldi.h:60
Array< OneD, NekDouble > m_imag_evl
Definition: DriverArnoldi.h:73
void CopyArnoldiArrayToField(Array< OneD, NekDouble > &array)
Copy Arnoldi storage to fields.
NekDouble m_evtol
Maxmum number of iterations.
Definition: DriverArnoldi.h:61
int m_nfields
interval to dump information if required.
Definition: DriverArnoldi.h:67
SOLVER_UTILS_EXPORT void ArnoldiSummary(std::ostream &out)
Array< OneD, NekDouble > m_real_evl
Operator in solve call is negated.
Definition: DriverArnoldi.h:72
void WriteEvs(std::ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble, bool DumpInverse=true)
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:89
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:86
enum EvolutionOperatorType m_EvolutionOperator
Evolution Operator.
Definition: Driver.h:104
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:98
int m_nequ
number of equations
Definition: Driver.h:101
virtual void v_InitObject(std::ostream &out=std::cout)
Virtual function for initialisation implementation.
void EV_update(Array< OneD, NekDouble > &src, Array< OneD, NekDouble > &tgt)
Generates a new vector in the sequence by applying the linear operator.
void EV_post(Array< OneD, Array< OneD, NekDouble > > &Tseq, Array< OneD, Array< OneD, NekDouble > > &Kseq, const int ntot, const int kdim, const int nvec, Array< OneD, NekDouble > &zvec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi, const int icon)
int EV_test(const int itrn, const int kdim, Array< OneD, NekDouble > &zvec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi, const NekDouble resnorm, int nvec, std::ofstream &evlout, NekDouble &resid0)
Tests for convergence of eigenvalues of H.
void EV_small(Array< OneD, Array< OneD, NekDouble > > &Kseq, Array< OneD, Array< OneD, NekDouble > > &Kseqcopy, const int ntot, const Array< OneD, NekDouble > &alpha, const int kdim, Array< OneD, NekDouble > &zvec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi, NekDouble &resnorm)
Generates the upper Hessenberg matrix H and computes its eigenvalues.
void EV_sort(Array< OneD, NekDouble > &evec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi, Array< OneD, NekDouble > &test, const int dim)
Sorts a sequence of eigenvectors/eigenvalues by magnitude.
void EV_big(Array< OneD, Array< OneD, NekDouble > > &bvecs, Array< OneD, Array< OneD, NekDouble > > &evecs, const int ntot, const int kdim, const int nvec, Array< OneD, NekDouble > &zvec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi)
virtual void v_Execute(std::ostream &out=std::cout)
Virtual function for solve implementation.
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition: Blas.hpp:197
std::shared_ptr< SessionReader > SessionReaderSharedPtr
DriverFactory & GetDriverFactory()
Definition: Driver.cpp:65
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:192
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:565
T Dot(int n, const T *w, const T *x)
vvtvp (vector times vector times vector): z = w*x*y
Definition: Vmath.cpp:1038
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:225
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436
void FillWhiteNoise(int n, const T eps, T *x, const int incx, int outseed)
Fills a vector with white noise.
Definition: Vmath.cpp:142
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267