Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | Static Private Attributes | Friends | List of all members
Nektar::SolverUtils::DriverModifiedArnoldi Class Reference

#include <DriverModifiedArnoldi.h>

Inheritance diagram for Nektar::SolverUtils::DriverModifiedArnoldi:
Inheritance graph
[legend]
Collaboration diagram for Nektar::SolverUtils::DriverModifiedArnoldi:
Collaboration graph
[legend]

Static Public Member Functions

static DriverSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of the class. More...
 

Protected Member Functions

 DriverModifiedArnoldi (const LibUtilities::SessionReaderSharedPtr pSession)
 Constructor. More...
 
virtual ~DriverModifiedArnoldi ()
 Destructor. More...
 
virtual void v_InitObject (ostream &out=cout)
 Virtual function for initialisation implementation. More...
 
virtual void v_Execute (ostream &out=cout)
 Virtual function for solve implementation. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::DriverArnoldi
 DriverArnoldi (const LibUtilities::SessionReaderSharedPtr pSession)
 Constructor. More...
 
virtual ~DriverArnoldi ()
 Destructor. More...
 
void CopyArnoldiArrayToField (Array< OneD, NekDouble > &array)
 Copy Arnoldi storage to fields. More...
 
void CopyFieldToArnoldiArray (Array< OneD, NekDouble > &array)
 Copy fields to Arnoldi storage. More...
 
void CopyFwdToAdj ()
 Copy the forward field to the adjoint system in transient growth calculations. More...
 
void WriteFld (std::string file, std::vector< Array< OneD, NekDouble > > coeffs)
 
void WriteFld (std::string file, Array< OneD, NekDouble > coeffs)
 
void WriteEvs (ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble)
 
virtual Array< OneD, NekDoublev_GetRealEvl (void)
 
virtual Array< OneD, NekDoublev_GetImagEvl (void)
 
- Protected Member Functions inherited from Nektar::SolverUtils::Driver
 Driver (const LibUtilities::SessionReaderSharedPtr pSession)
 Initialises EquationSystem class members. More...
 

Private Member Functions

void EV_update (Array< OneD, NekDouble > &src, Array< OneD, NekDouble > &tgt)
 Generates a new vector in the sequence by applying the linear operator. More...
 
void EV_small (Array< OneD, Array< OneD, NekDouble > > &Kseq, 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. More...
 
int EV_test (const int itrn, const int kdim, Array< OneD, NekDouble > &zvec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi, const NekDouble resnorm, const int nvec, ofstream &evlout, NekDouble &resid0)
 Tests for convergence of eigenvalues of H. More...
 
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. More...
 
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)
 
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)
 

Static Private Attributes

static std::string driverLookupId
 

Friends

class MemoryManager< DriverModifiedArnoldi >
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::SolverUtils::DriverArnoldi
SOLVER_UTILS_EXPORT void ArnoldiSummary (std::ostream &out)
 
- Public Member Functions inherited from Nektar::SolverUtils::Driver
virtual ~Driver ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void InitObject (ostream &out=cout)
 Initialise Object. More...
 
SOLVER_UTILS_EXPORT void Execute (ostream &out=cout)
 Execute driver. More...
 
SOLVER_UTILS_EXPORT Array< OneD, EquationSystemSharedPtrGetEqu ()
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleGetRealEvl (void)
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleGetImagEvl (void)
 
- Protected Attributes inherited from Nektar::SolverUtils::DriverArnoldi
int m_kdim
 
int m_nvec
 Dimension of Krylov subspace. More...
 
int m_nits
 Number of vectors to test. More...
 
NekDouble m_evtol
 Maxmum number of iterations. More...
 
NekDouble m_period
 Tolerance of iteratiosn. More...
 
bool m_timeSteppingAlgorithm
 Period of time stepping algorithm. More...
 
int m_infosteps
 underlying operator is time stepping More...
 
int m_nfields
 interval to dump information if required. More...
 
NekDouble m_realShift
 
NekDouble m_imagShift
 
Array< OneD, NekDoublem_real_evl
 
Array< OneD, NekDoublem_imag_evl
 
- Protected Attributes inherited from Nektar::SolverUtils::Driver
LibUtilities::CommSharedPtr m_comm
 Communication object. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 Session reader object. More...
 
LibUtilities::SessionReaderSharedPtr session_LinNS
 I the Coupling between SFD and arnoldi. More...
 
Array< OneD, EquationSystemSharedPtrm_equ
 Equation system to solve. More...
 
int m_nequ
 number of equations More...
 
enum EvolutionOperatorType m_EvolutionOperator
 Evolution Operator. More...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::Driver
static std::string evolutionOperatorLookupIds []
 
static std::string evolutionOperatorDef
 
static std::string driverDefault
 

Detailed Description

Definition at line 47 of file DriverModifiedArnoldi.h.

Constructor & Destructor Documentation

Nektar::SolverUtils::DriverModifiedArnoldi::DriverModifiedArnoldi ( const LibUtilities::SessionReaderSharedPtr  pSession)
protected

Constructor.

Definition at line 56 of file DriverModifiedArnoldi.cpp.

58  : DriverArnoldi(pSession)
59  {
60  }
DriverArnoldi(const LibUtilities::SessionReaderSharedPtr pSession)
Constructor.
Nektar::SolverUtils::DriverModifiedArnoldi::~DriverModifiedArnoldi ( )
protectedvirtual

Destructor.

Definition at line 66 of file DriverModifiedArnoldi.cpp.

67  {
68  }

Member Function Documentation

static DriverSharedPtr Nektar::SolverUtils::DriverModifiedArnoldi::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Definition at line 53 of file DriverModifiedArnoldi.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

55  {
57  ::AllocateSharedPtr(pSession);
58  p->InitObject();
59  return p;
60  }
boost::shared_ptr< Driver > DriverSharedPtr
A shared pointer to a Driver object.
Definition: Driver.h:52
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
void Nektar::SolverUtils::DriverModifiedArnoldi::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 
)
private

Definition at line 527 of file DriverModifiedArnoldi.cpp.

References Vmath::Ddot(), Nektar::SolverUtils::Driver::m_comm, Nektar::LibUtilities::ReduceSum, Vmath::Smul(), Vmath::Svtvp(), and Vmath::Zero().

Referenced by EV_post().

536  {
537  NekDouble wgt, norm;
538 
539  // Generate big eigenvectors
540  for (int j = 0; j < nvec; ++j)
541  {
542  Vmath::Zero(ntot, evecs[j], 1);
543  for (int i = 0; i < kdim; ++i)
544  {
545  wgt = zvec[i + j*kdim];
546  Vmath::Svtvp(ntot, wgt, bvecs[i], 1, evecs[j], 1, evecs[j], 1);
547  }
548  }
549 
550  // Normalise the big eigenvectors
551  for (int i = 0; i < nvec; ++i)
552  {
553  if (wi[i] == 0.0) // Real mode
554  {
555  norm = Blas::Ddot(ntot, &evecs[i][0], 1, &evecs[i][0], 1);
556  m_comm->AllReduce(norm, Nektar::LibUtilities::ReduceSum);
557  norm = std::sqrt(norm);
558  Vmath::Smul(ntot, 1.0/norm, evecs[i], 1, evecs[i], 1);
559  }
560  else
561  {
562  norm = Blas::Ddot(ntot, &evecs[i][0], 1, &evecs[i][0], 1);
563  norm += Blas::Ddot(ntot, &evecs[i+1][0], 1, &evecs[i+1][0], 1);
564  m_comm->AllReduce(norm, Nektar::LibUtilities::ReduceSum);
565  norm = std::sqrt(norm);
566 
567  Vmath::Smul(ntot, 1.0/norm, evecs[i], 1, evecs[i], 1);
568  Vmath::Smul(ntot, 1.0/norm, evecs[i+1], 1, evecs[i+1], 1);
569 
570  i++;
571  }
572  }
573  }
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:471
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
double NekDouble
T Ddot(int n, const Array< OneD, const T > &w, const int incw, const Array< OneD, const T > &x, const int incx, const Array< OneD, const int > &y, const int incy)
Definition: VmathArray.hpp:426
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:85
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
void Nektar::SolverUtils::DriverModifiedArnoldi::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 
)
private

Definition at line 479 of file DriverModifiedArnoldi.cpp.

References ASSERTL0, EV_big(), Nektar::SolverUtils::Driver::m_equ, Nektar::SolverUtils::Driver::m_session, and Nektar::SolverUtils::DriverArnoldi::WriteFld().

Referenced by v_Execute().

489  {
490  if (icon == 0)
491  {
492  // Not converged, write final Krylov vector
493  ASSERTL0(false, "Convergence was not achieved within the "
494  "prescribed number of iterations.");
495  }
496  else if (icon < 0)
497  {
498  // Minimum residual reached
499  ASSERTL0(false, "Minimum residual reached.");
500  }
501  else if (icon == nvec)
502  {
503  // Converged, write out eigenvectors
504  EV_big(Tseq, Kseq, ntot, kdim, icon, zvec, wr, wi);
506  = m_equ[0]->UpdateFields();
507 
508  for (int j = 0; j < icon; ++j)
509  {
510  std::string file = m_session->GetSessionName() + "_eig_"
511  + boost::lexical_cast<std::string>(j);
512 
513  WriteFld(file,Kseq[j]);
514  }
515  }
516  else
517  {
518  // Not recognised
519  ASSERTL0(false, "Unrecognised value.");
520  }
521  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
void WriteFld(std::string file, std::vector< Array< OneD, NekDouble > > coeffs)
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:94
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)
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:88
void Nektar::SolverUtils::DriverModifiedArnoldi::EV_small ( Array< OneD, Array< OneD, NekDouble > > &  Kseq,
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 
)
private

Generates the upper Hessenberg matrix H and computes its eigenvalues.

Definition at line 317 of file DriverModifiedArnoldi.cpp.

References ASSERTL0, Vmath::Ddot(), Vmath::Dot(), Nektar::SolverUtils::Driver::m_comm, Nektar::LibUtilities::ReduceSum, Vmath::Smul(), and Vmath::Svtvp().

Referenced by v_Execute().

326  {
327  int kdimp = kdim + 1;
328  int lwork = 10*kdim;
329  int ier;
330  Array<OneD, NekDouble> R(kdimp * kdimp, 0.0);
331  Array<OneD, NekDouble> H(kdimp * kdim, 0.0);
332  Array<OneD, NekDouble> rwork(lwork, 0.0);
333 
334  // Modified G-S orthonormalisation
335  for (int i = 0; i < kdimp; ++i)
336  {
337  NekDouble gsc = Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[i][0], 1);
338  m_comm->AllReduce(gsc, Nektar::LibUtilities::ReduceSum);
339  gsc = std::sqrt(gsc);
340  ASSERTL0(gsc != 0.0, "Vectors are linearly independent.");
341 
342  R[i*kdimp+i] = gsc;
343  Vmath::Smul(ntot, 1.0/gsc, Kseq[i], 1, Kseq[i], 1);
344 
345  for (int j = i + 1; j < kdimp; ++j)
346  {
347  gsc = Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[j][0], 1);
348  m_comm->AllReduce(gsc, Nektar::LibUtilities::ReduceSum);
349  Vmath::Svtvp(ntot, -gsc, Kseq[i], 1, Kseq[j], 1, Kseq[j], 1);
350  R[j*kdimp+i] = gsc;
351  }
352  }
353 
354  // Compute H matrix
355  for (int i = 0; i < kdim; ++i)
356  {
357  for (int j = 0; j < kdim; ++j)
358  {
359  H[j*kdim+i] = alpha[j+1] * R[(j+1)*kdimp+i]
360  - Vmath::Dot(j, &H[0] + i, kdim, &R[0] + j*kdimp, 1);
361  H[j*kdim+i] /= R[j*kdimp+j];
362  }
363  }
364 
365  H[(kdim-1)*kdim+kdim] = alpha[kdim]
366  * std::fabs(R[kdim*kdimp+kdim] / R[(kdim-1)*kdimp + kdim-1]);
367 
368  Lapack::dgeev_('N', 'V', kdim, &H[0], kdim, &wr[0], &wi[0], 0, 1,
369  &zvec[0], kdim, &rwork[0], lwork, ier);
370 
371  ASSERTL0(!ier, "Error with dgeev");
372 
373  resnorm = H[(kdim-1)*kdim + kdim];
374  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
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:471
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
double NekDouble
T Ddot(int n, const Array< OneD, const T > &w, const int incw, const Array< OneD, const T > &x, const int incx, const Array< OneD, const int > &y, const int incy)
Definition: VmathArray.hpp:426
T Dot(int n, const T *w, const T *x)
vvtvp (vector times vector times vector): z = w*x*y
Definition: Vmath.cpp:869
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:85
void Nektar::SolverUtils::DriverModifiedArnoldi::EV_sort ( Array< OneD, NekDouble > &  evec,
Array< OneD, NekDouble > &  wr,
Array< OneD, NekDouble > &  wi,
Array< OneD, NekDouble > &  test,
const int  dim 
)
private

Sorts a sequence of eigenvectors/eigenvalues by magnitude.

Definition at line 444 of file DriverModifiedArnoldi.cpp.

References Vmath::Vcopy().

Referenced by EV_test().

450  {
451  Array<OneD, NekDouble> z_tmp(dim,0.0);
452  NekDouble wr_tmp, wi_tmp, te_tmp;
453  for (int j = 1; j < dim; ++j)
454  {
455  wr_tmp = wr[j];
456  wi_tmp = wi[j];
457  te_tmp = test[j];
458  Vmath::Vcopy(dim, &evec[0] + j*dim, 1, &z_tmp[0], 1);
459  int i = j - 1;
460  while (i >= 0 && test[i] > te_tmp)
461  {
462  wr[i+1] = wr[i];
463  wi[i+1] = wi[i];
464  test[i+1] = test[i];
465  Vmath::Vcopy(dim, &evec[0] + i*dim, 1, &evec[0] + (i+1)*dim, 1);
466  i--;
467  }
468  wr[i+1] = wr_tmp;
469  wi[i+1] = wi_tmp;
470  test[i+1] = te_tmp;
471  Vmath::Vcopy(dim, &z_tmp[0], 1, &evec[0] + (i+1)*dim, 1);
472  }
473  }
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1016
int Nektar::SolverUtils::DriverModifiedArnoldi::EV_test ( const int  itrn,
const int  kdim,
Array< OneD, NekDouble > &  zvec,
Array< OneD, NekDouble > &  wr,
Array< OneD, NekDouble > &  wi,
const NekDouble  resnorm,
const int  nvec,
ofstream &  evlout,
NekDouble resid0 
)
private

Tests for convergence of eigenvalues of H.

Definition at line 380 of file DriverModifiedArnoldi.cpp.

References Vmath::Dot(), EV_sort(), Nektar::SolverUtils::Driver::m_comm, Nektar::SolverUtils::DriverArnoldi::m_evtol, Nektar::SolverUtils::DriverArnoldi::m_timeSteppingAlgorithm, and Nektar::SolverUtils::DriverArnoldi::WriteEvs().

Referenced by v_Execute().

390  {
391  int idone = 0;
392  // NekDouble period = 0.1;
393 
394  Array<OneD, NekDouble> resid(kdim);
395  for (int i = 0; i < kdim; ++i)
396  {
397  NekDouble tmp = std::sqrt(Vmath::Dot(kdim, &zvec[0] + i*kdim, 1,
398  &zvec[0] + i*kdim, 1));
399  resid[i] = resnorm * std::fabs(zvec[kdim - 1 + i*kdim]) / tmp;
400  if (wi[i] < 0.0)
401  {
402  resid[i-1] = resid[i] = hypot(resid[i-1], resid[i]);
403  }
404  }
405 
406  EV_sort(zvec, wr, wi, resid, kdim);
407 
408  if (resid[nvec-1] < m_evtol)
409  {
410  idone = nvec;
411  }
412 
413  if (m_comm->GetRank() == 0)
414  {
415  evlout << "-- Iteration = " << itrn << ", H(k+1, k) = "
416  << resnorm << endl;
417  evlout.precision(4);
418  evlout.setf(ios::scientific, ios::floatfield);
420  {
421  evlout << "EV Magnitude Angle Growth "
422  << "Frequency Residual" << endl;
423  }
424  else
425  {
426  evlout << "EV Real Imaginary inverse real "
427  << "inverse imag Residual" << endl;
428  }
429 
430  for (int i = 0; i < kdim; i++)
431  {
432  WriteEvs(evlout,i,wr[i],wi[i],resid[i]);
433  }
434  }
435 
436  resid0 = resid[nvec-1];
437  return idone;
438  }
NekDouble m_evtol
Maxmum number of iterations.
Definition: DriverArnoldi.h:57
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 WriteEvs(ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble)
double NekDouble
bool m_timeSteppingAlgorithm
Period of time stepping algorithm.
Definition: DriverArnoldi.h:59
T Dot(int n, const T *w, const T *x)
vvtvp (vector times vector times vector): z = w*x*y
Definition: Vmath.cpp:869
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:85
void Nektar::SolverUtils::DriverModifiedArnoldi::EV_update ( Array< OneD, NekDouble > &  src,
Array< OneD, NekDouble > &  tgt 
)
private

Generates a new vector in the sequence by applying the linear operator.

Definition at line 287 of file DriverModifiedArnoldi.cpp.

References Nektar::SolverUtils::DriverArnoldi::CopyArnoldiArrayToField(), Nektar::SolverUtils::DriverArnoldi::CopyFieldToArnoldiArray(), Nektar::SolverUtils::DriverArnoldi::CopyFwdToAdj(), Nektar::SolverUtils::eTransientGrowth, Nektar::SolverUtils::Driver::m_equ, and Nektar::SolverUtils::Driver::m_EvolutionOperator.

Referenced by v_Execute().

290  {
291  // Copy starting vector into first sequence element.
293  m_equ[0]->TransCoeffToPhys();
294 
295  m_equ[0]->DoSolve();
296 
298  {
300  fields = m_equ[0]->UpdateFields();
301 
302  //start Adjoint with latest fields of direct
303  CopyFwdToAdj();
304  m_equ[1]->TransCoeffToPhys();
305 
306  m_equ[1]->DoSolve();
307  }
308 
309  // Copy starting vector into first sequence element.
311  }
void CopyFwdToAdj()
Copy the forward field to the adjoint system in transient growth calculations.
void CopyArnoldiArrayToField(Array< OneD, NekDouble > &array)
Copy Arnoldi storage to fields.
void CopyFieldToArnoldiArray(Array< OneD, NekDouble > &array)
Copy fields to Arnoldi storage.
enum EvolutionOperatorType m_EvolutionOperator
Evolution Operator.
Definition: Driver.h:100
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:94
void Nektar::SolverUtils::DriverModifiedArnoldi::v_Execute ( ostream &  out = cout)
protectedvirtual

Virtual function for solve implementation.

Implements Nektar::SolverUtils::Driver.

Reimplemented in Nektar::SolverUtils::DriverSteadyState.

Definition at line 99 of file DriverModifiedArnoldi.cpp.

References Nektar::SolverUtils::DriverArnoldi::CopyFieldToArnoldiArray(), Vmath::Ddot(), EV_post(), EV_small(), EV_test(), EV_update(), Vmath::FillWhiteNoise(), Nektar::SolverUtils::Driver::m_comm, Nektar::SolverUtils::Driver::m_equ, Nektar::SolverUtils::DriverArnoldi::m_imag_evl, Nektar::SolverUtils::DriverArnoldi::m_kdim, Nektar::SolverUtils::DriverArnoldi::m_nfields, Nektar::SolverUtils::DriverArnoldi::m_nits, Nektar::SolverUtils::DriverArnoldi::m_nvec, Nektar::SolverUtils::DriverArnoldi::m_real_evl, Nektar::SolverUtils::Driver::m_session, Nektar::LibUtilities::ReduceSum, Vmath::Smul(), and Vmath::Vcopy().

Referenced by Nektar::SolverUtils::DriverSteadyState::v_Execute().

100  {
101  int i = 0;
102  int j = 0;
103  int nq = m_equ[0]->UpdateFields()[0]->GetNcoeffs();
104  int ntot = m_nfields*nq;
105  int converged = 0;
106  NekDouble resnorm = 0.0;
107  ofstream evlout;
108  std::string evlFile = m_session->GetSessionName() + ".evl";
109 
110  if (m_comm->GetRank() == 0)
111  {
112  evlout.open(evlFile.c_str());
113  }
114 
115  // Allocate memory
120 
125  for (i = 0; i < m_kdim + 1; ++i)
126  {
127  Kseq[i] = Array<OneD, NekDouble>(ntot, 0.0);
128  Tseq[i] = Array<OneD, NekDouble>(ntot, 0.0);
129  }
130 
131  // Copy starting vector into second sequence element (temporary).
132  if(m_session->DefinesFunction("InitialConditions"))
133  {
134  if (m_comm->GetRank() == 0)
135  {
136  out << "\tInital vector : specified in input file " << endl;
137  }
138  m_equ[0]->SetInitialConditions(0.0,false);
139 
140  CopyFieldToArnoldiArray(Kseq[1]);
141  }
142  else
143  {
144  if (m_comm->GetRank() == 0)
145  {
146  out << "\tInital vector : random " << endl;
147  }
148 
149  NekDouble eps=1;
150  Vmath::FillWhiteNoise(ntot, eps , &Kseq[1][0], 1);
151  }
152 
153  // Perform one iteration to enforce boundary conditions.
154  // Set this as the initial value in the sequence.
155  EV_update(Kseq[1], Kseq[0]);
156  if (m_comm->GetRank() == 0)
157  {
158  out << "Iteration: " << 0 << endl;
159  }
160 
161  // Normalise first vector in sequence
162  alpha[0] = Blas::Ddot(ntot, &Kseq[0][0], 1, &Kseq[0][0], 1);
163  m_comm->AllReduce(alpha[0], Nektar::LibUtilities::ReduceSum);
164  alpha[0] = std::sqrt(alpha[0]);
165  Vmath::Smul(ntot, 1.0/alpha[0], Kseq[0], 1, Kseq[0], 1);
166 
167  // Fill initial krylov sequence
168  NekDouble resid0;
169  for (i = 1; !converged && i <= m_kdim; ++i)
170  {
171  // Compute next vector
172  EV_update(Kseq[i-1], Kseq[i]);
173 
174  // Normalise
175  alpha[i] = Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[i][0], 1);
176  m_comm->AllReduce(alpha[i], Nektar::LibUtilities::ReduceSum);
177  alpha[i] = std::sqrt(alpha[i]);
178 
179  //alpha[i] = std::sqrt(alpha[i]);
180  Vmath::Smul(ntot, 1.0/alpha[i], Kseq[i], 1, Kseq[i], 1);
181 
182  // Copy Krylov sequence into temporary storage
183  for (int k = 0; k < i + 1; ++k)
184  {
185  Vmath::Vcopy(ntot, Kseq[k], 1, Tseq[k], 1);
186  }
187 
188  // Generate Hessenberg matrix and compute eigenvalues of it.
189  EV_small(Tseq, ntot, alpha, i, zvec, wr, wi, resnorm);
190 
191  // Test for convergence.
192  converged = EV_test(i, i, zvec, wr, wi, resnorm,
193  std::min(i, m_nvec), evlout, resid0);
194  converged = max (converged, 0);
195 
196  if (m_comm->GetRank() == 0)
197  {
198  out << "Iteration: " << i << " (residual : " << resid0
199  << ")" <<endl;
200  }
201  }
202 
203  // Continue with full sequence
204  if (!converged)
205  {
206  for (i = m_kdim + 1; !converged && i <= m_nits; ++i)
207  {
208  // Shift all the vectors in the sequence.
209  // First vector is removed.
210  for (int j = 1; j <= m_kdim; ++j)
211  {
212  alpha[j-1] = alpha[j];
213  Vmath::Vcopy(ntot, Kseq[j], 1, Kseq[j-1], 1);
214  }
215 
216  // Compute next vector
217  EV_update(Kseq[m_kdim - 1], Kseq[m_kdim]);
218 
219  // Compute new scale factor
220  alpha[m_kdim] = Blas::Ddot(ntot, &Kseq[m_kdim][0], 1,
221  &Kseq[m_kdim][0], 1);
222  m_comm->AllReduce(alpha[m_kdim],
224  alpha[m_kdim] = std::sqrt(alpha[m_kdim]);
225  Vmath::Smul(ntot, 1.0/alpha[m_kdim], Kseq[m_kdim], 1,
226  Kseq[m_kdim], 1);
227 
228  // Copy Krylov sequence into temporary storage
229  for (int k = 0; k < m_kdim + 1; ++k)
230  {
231  Vmath::Vcopy(ntot, Kseq[k], 1, Tseq[k], 1);
232  }
233 
234  // Generate Hessenberg matrix and compute eigenvalues of it
235  EV_small(Tseq, ntot, alpha, m_kdim, zvec, wr, wi, resnorm);
236 
237  // Test for convergence.
238  converged = EV_test(i, m_kdim, zvec, wr, wi, resnorm,
239  m_nvec, evlout, resid0);
240 
241  if (m_comm->GetRank() == 0)
242  {
243  out << "Iteration: " << i << " (residual : "
244  << resid0 << ")" <<endl;
245  }
246  }
247  }
248 
249  m_equ[0]->Output();
250 
251  // Evaluate and output computation time and solution accuracy.
252  // The specific format of the error output is essential for the
253  // regression tests to work.
254  // Evaluate L2 Error
255  for(j = 0; j < m_equ[0]->GetNvariables(); ++j)
256  {
257  NekDouble vL2Error = m_equ[0]->L2Error(j,false);
258  NekDouble vLinfError = m_equ[0]->LinfError(j);
259  if (m_comm->GetRank() == 0)
260  {
261  out << "L 2 error (variable " << m_equ[0]->GetVariable(j)
262  << ") : " << vL2Error << endl;
263  out << "L inf error (variable " << m_equ[0]->GetVariable(j)
264  << ") : " << vLinfError << endl;
265  }
266  }
267 
268  // Process eigenvectors and write out.
269  EV_post(Tseq, Kseq, ntot, min(--i, m_kdim), m_nvec, zvec, wr, wi,
270  converged);
271 
272  // store eigenvalues so they can be access from driver class
273  m_real_evl = wr;
274  m_imag_evl = wi;
275 
276  // Close the runtime info file.
277  if (m_comm->GetRank() == 0)
278  {
279  evlout.close();
280  }
281  }
Array< OneD, NekDouble > m_imag_evl
Definition: DriverArnoldi.h:68
void EV_small(Array< OneD, Array< OneD, NekDouble > > &Kseq, 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 CopyFieldToArnoldiArray(Array< OneD, NekDouble > &array)
Copy fields to Arnoldi storage.
void EV_update(Array< OneD, NekDouble > &src, Array< OneD, NekDouble > &tgt)
Generates a new vector in the sequence by applying the linear operator.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
int m_nvec
Dimension of Krylov subspace.
Definition: DriverArnoldi.h:55
Array< OneD, NekDouble > m_real_evl
Definition: DriverArnoldi.h:67
double NekDouble
T Ddot(int n, const Array< OneD, const T > &w, const int incw, const Array< OneD, const T > &x, const int incx, const Array< OneD, const int > &y, const int incy)
Definition: VmathArray.hpp:426
int m_nfields
interval to dump information if required.
Definition: DriverArnoldi.h:63
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:94
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)
void FillWhiteNoise(int n, const T eps, T *x, const int incx, int outseed)
Fills a vector with white noise.
Definition: Vmath.cpp:138
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:85
int EV_test(const int itrn, const int kdim, Array< OneD, NekDouble > &zvec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi, const NekDouble resnorm, const int nvec, ofstream &evlout, NekDouble &resid0)
Tests for convergence of eigenvalues of H.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1016
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:88
int m_nits
Number of vectors to test.
Definition: DriverArnoldi.h:56
void Nektar::SolverUtils::DriverModifiedArnoldi::v_InitObject ( ostream &  out = cout)
protectedvirtual

Virtual function for initialisation implementation.

Reimplemented from Nektar::SolverUtils::DriverArnoldi.

Reimplemented in Nektar::SolverUtils::DriverSteadyState.

Definition at line 74 of file DriverModifiedArnoldi.cpp.

References Nektar::SolverUtils::DriverArnoldi::ArnoldiSummary(), Nektar::SolverUtils::Driver::m_comm, Nektar::SolverUtils::Driver::m_equ, Nektar::SolverUtils::Driver::m_nequ, and Nektar::SolverUtils::DriverArnoldi::v_InitObject().

Referenced by Nektar::SolverUtils::DriverSteadyState::v_InitObject().

75  {
77 
78  m_equ[0]->PrintSummary(out);
79 
80  // Print session parameters
81  if (m_comm->GetRank() == 0)
82  {
83  out << "\tArnoldi solver type : Modified Arnoldi" << endl;
84  }
85 
87 
88  m_equ[m_nequ - 1]->DoInitialise();
89 
90  //FwdTrans Initial conditions to be in Coefficient Space
91  m_equ[m_nequ-1] ->TransPhysToCoeff();
92 
93  }
SOLVER_UTILS_EXPORT void ArnoldiSummary(std::ostream &out)
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:94
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:85
virtual void v_InitObject(ostream &out=cout)
int m_nequ
number of equations
Definition: Driver.h:97

Friends And Related Function Documentation

friend class MemoryManager< DriverModifiedArnoldi >
friend

Definition at line 50 of file DriverModifiedArnoldi.h.

Member Data Documentation

string Nektar::SolverUtils::DriverModifiedArnoldi::className
static
Initial value:

Name of the class.

Definition at line 63 of file DriverModifiedArnoldi.h.

string Nektar::SolverUtils::DriverModifiedArnoldi::driverLookupId
staticprivate
Initial value:
=
"ModifiedArnoldi",0)

Definition at line 132 of file DriverModifiedArnoldi.h.