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:
[legend]

Static Public Member Functions

static DriverSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 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, const SpatialDomains::MeshGraphSharedPtr pGraph)
 Constructor. More...
 
virtual ~DriverModifiedArnoldi ()
 Destructor. More...
 
virtual void v_InitObject (std::ostream &out=std::cout)
 Virtual function for initialisation implementation. More...
 
virtual void v_Execute (std::ostream &out=std::cout)
 Virtual function for solve implementation. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::DriverArnoldi
 DriverArnoldi (const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
 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)
 Write coefficients to file. More...
 
void WriteFld (std::string file, Array< OneD, NekDouble > coeffs)
 
void WriteEvs (std::ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble, bool DumpInverse=true)
 
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, const SpatialDomains::MeshGraphSharedPtr pGraph)
 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, std::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 (std::ostream &out=std::cout)
 Initialise Object. More...
 
SOLVER_UTILS_EXPORT void Execute (std::ostream &out=std::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
 
int m_negatedOp
 
Array< OneD, NekDoublem_real_evl
 Operator in solve call is negated. More...
 
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...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 MeshGraph object. 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 46 of file DriverModifiedArnoldi.h.

Constructor & Destructor Documentation

◆ DriverModifiedArnoldi()

Nektar::SolverUtils::DriverModifiedArnoldi::DriverModifiedArnoldi ( const LibUtilities::SessionReaderSharedPtr  pSession,
const SpatialDomains::MeshGraphSharedPtr  pGraph 
)
protected

Constructor.

Definition at line 59 of file DriverModifiedArnoldi.cpp.

62  : DriverArnoldi(pSession, pGraph)
63 {
64 }
DriverArnoldi(const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
Constructor.

◆ ~DriverModifiedArnoldi()

Nektar::SolverUtils::DriverModifiedArnoldi::~DriverModifiedArnoldi ( )
protectedvirtual

Destructor.

Definition at line 70 of file DriverModifiedArnoldi.cpp.

71 {
72 }

Member Function Documentation

◆ create()

static DriverSharedPtr Nektar::SolverUtils::DriverModifiedArnoldi::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
inlinestatic

Creates an instance of this class.

Definition at line 52 of file DriverModifiedArnoldi.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

55  {
57  ::AllocateSharedPtr(pSession, pGraph);
58  p->InitObject();
59  return p;
60  }
std::shared_ptr< Driver > DriverSharedPtr
A shared pointer to a Driver object.
Definition: Driver.h:51
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

◆ EV_big()

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

Compute estimates of the eigenvalues/eigenvectors of the original system.

Definition at line 555 of file DriverModifiedArnoldi.cpp.

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

Referenced by EV_post().

564 {
565  boost::ignore_unused(wr);
566 
567  NekDouble wgt, norm;
568 
569  // Generate big eigenvectors
570  for (int j = 0; j < nvec; ++j)
571  {
572  Vmath::Zero(ntot, evecs[j], 1);
573  for (int i = 0; i < kdim; ++i)
574  {
575  wgt = zvec[i + j*kdim];
576  Vmath::Svtvp(ntot, wgt, bvecs[i], 1, evecs[j], 1, evecs[j], 1);
577  }
578  }
579 
580  // Normalise the big eigenvectors
581  for (int i = 0; i < nvec; ++i)
582  {
583  if (wi[i] == 0.0) // Real mode
584  {
585  norm = Blas::Ddot(ntot, &evecs[i][0], 1, &evecs[i][0], 1);
586  m_comm->AllReduce(norm, Nektar::LibUtilities::ReduceSum);
587  norm = std::sqrt(norm);
588  Vmath::Smul(ntot, 1.0/norm, evecs[i], 1, evecs[i], 1);
589  }
590  else
591  {
592  norm = Blas::Ddot(ntot, &evecs[i][0], 1, &evecs[i][0], 1);
593  norm += Blas::Ddot(ntot, &evecs[i+1][0], 1, &evecs[i+1][0], 1);
594  m_comm->AllReduce(norm, Nektar::LibUtilities::ReduceSum);
595  norm = std::sqrt(norm);
596 
597  Vmath::Smul(ntot, 1.0/norm, evecs[i], 1, evecs[i], 1);
598  Vmath::Smul(ntot, 1.0/norm, evecs[i+1], 1, evecs[i+1], 1);
599 
600  i++;
601  }
602  }
603 }
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:488
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:216
double NekDouble
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:140
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:86
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376

◆ EV_post()

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

Post-process the Ritz eigenvalues/eigenvectors of the matrix H, to compute estimations of the leading eigenvalues and eigenvectors of the original matrix.

Definition at line 502 of file DriverModifiedArnoldi.cpp.

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

Referenced by v_Execute().

512 {
513  if (icon == 0)
514  {
515  // Not converged, write final Krylov vector
516  ASSERTL0(false, "Convergence was not achieved within the "
517  "prescribed number of iterations.");
518  }
519  else if (icon < 0)
520  {
521  // Minimum residual reached
522  ASSERTL0(false, "Minimum residual reached.");
523  }
524  else if (icon == nvec)
525  {
526  // Converged, write out eigenvectors
527  EV_big(Tseq, Kseq, ntot, kdim, icon, zvec, wr, wi);
528  Array<OneD, MultiRegions::ExpListSharedPtr> fields
529  = m_equ[0]->UpdateFields();
530 
531  for (int j = 0; j < icon; ++j)
532  {
533  std::string file = m_session->GetSessionName() + "_eig_"
534  + boost::lexical_cast<std::string>(j)
535  + ".fld";
536 
537  if (m_comm->GetRank() == 0)
538  {
539  WriteEvs(cout, j, wr[j], wi[j]);
540  }
541  WriteFld(file,Kseq[j]);
542  }
543  }
544  else
545  {
546  // Not recognised
547  ASSERTL0(false, "Unrecognised value.");
548  }
549 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
void WriteFld(std::string file, std::vector< Array< OneD, NekDouble > > coeffs)
Write coefficients to file.
void WriteEvs(std::ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble, bool DumpInverse=true)
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:98
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::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:86
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:89

◆ EV_small()

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.

Computes the Ritz eigenvalues and eigenvectors of the Hessenberg matrix constructed from the Krylov sequence.

Definition at line 338 of file DriverModifiedArnoldi.cpp.

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

Referenced by v_Execute().

347 {
348  int kdimp = kdim + 1;
349  int lwork = 10*kdim;
350  int ier;
351  Array<OneD, NekDouble> R(kdimp * kdimp, 0.0);
352  Array<OneD, NekDouble> H(kdimp * kdim, 0.0);
353  Array<OneD, NekDouble> rwork(lwork, 0.0);
354 
355  // Modified G-S orthonormalisation
356  for (int i = 0; i < kdimp; ++i)
357  {
358  NekDouble gsc = Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[i][0], 1);
359  m_comm->AllReduce(gsc, Nektar::LibUtilities::ReduceSum);
360  gsc = std::sqrt(gsc);
361  ASSERTL0(gsc != 0.0, "Vectors are linearly independent.");
362 
363  R[i*kdimp+i] = gsc;
364  Vmath::Smul(ntot, 1.0/gsc, Kseq[i], 1, Kseq[i], 1);
365 
366  for (int j = i + 1; j < kdimp; ++j)
367  {
368  gsc = Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[j][0], 1);
369  m_comm->AllReduce(gsc, Nektar::LibUtilities::ReduceSum);
370  Vmath::Svtvp(ntot, -gsc, Kseq[i], 1, Kseq[j], 1, Kseq[j], 1);
371  R[j*kdimp+i] = gsc;
372  }
373  }
374 
375  // Compute H matrix
376  for (int i = 0; i < kdim; ++i)
377  {
378  for (int j = 0; j < kdim; ++j)
379  {
380  H[j*kdim+i] = alpha[j+1] * R[(j+1)*kdimp+i]
381  - Vmath::Dot(j, &H[0] + i, kdim, &R[0] + j*kdimp, 1);
382  H[j*kdim+i] /= R[j*kdimp+j];
383  }
384  }
385 
386  H[(kdim-1)*kdim+kdim] = alpha[kdim]
387  * std::fabs(R[kdim*kdimp+kdim] / R[(kdim-1)*kdimp + kdim-1]);
388 
389  Lapack::dgeev_('N', 'V', kdim, &H[0], kdim, &wr[0], &wi[0], 0, 1,
390  &zvec[0], kdim, &rwork[0], lwork, ier);
391 
392  ASSERTL0(!ier, "Error with dgeev");
393 
394  resnorm = H[(kdim-1)*kdim + kdim];
395 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
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:488
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:216
double NekDouble
T Dot(int n, const T *w, const T *x)
vvtvp (vector times vector times vector): z = w*x*y
Definition: Vmath.cpp:917
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:140
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:86

◆ EV_sort()

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.

Sorts the computed eigenvalues by smallest residual first.

Definition at line 465 of file DriverModifiedArnoldi.cpp.

References Vmath::Vcopy().

Referenced by EV_test().

471 {
472  Array<OneD, NekDouble> z_tmp(dim,0.0);
473  NekDouble wr_tmp, wi_tmp, te_tmp;
474  for (int j = 1; j < dim; ++j)
475  {
476  wr_tmp = wr[j];
477  wi_tmp = wi[j];
478  te_tmp = test[j];
479  Vmath::Vcopy(dim, &evec[0] + j*dim, 1, &z_tmp[0], 1);
480  int i = j - 1;
481  while (i >= 0 && test[i] > te_tmp)
482  {
483  wr[i+1] = wr[i];
484  wi[i+1] = wi[i];
485  test[i+1] = test[i];
486  Vmath::Vcopy(dim, &evec[0] + i*dim, 1, &evec[0] + (i+1)*dim, 1);
487  i--;
488  }
489  wr[i+1] = wr_tmp;
490  wi[i+1] = wi_tmp;
491  test[i+1] = te_tmp;
492  Vmath::Vcopy(dim, &z_tmp[0], 1, &evec[0] + (i+1)*dim, 1);
493  }
494 }
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ EV_test()

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,
std::ofstream &  evlout,
NekDouble resid0 
)
private

Tests for convergence of eigenvalues of H.

Check for convergence of the residuals of the eigenvalues of H.

Definition at line 401 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().

411 {
412  int idone = 0;
413  // NekDouble period = 0.1;
414 
415  Array<OneD, NekDouble> resid(kdim);
416  for (int i = 0; i < kdim; ++i)
417  {
418  NekDouble tmp = std::sqrt(Vmath::Dot(kdim, &zvec[0] + i*kdim, 1,
419  &zvec[0] + i*kdim, 1));
420  resid[i] = resnorm * std::fabs(zvec[kdim - 1 + i*kdim]) / tmp;
421  if (wi[i] < 0.0)
422  {
423  resid[i-1] = resid[i] = hypot(resid[i-1], resid[i]);
424  }
425  }
426 
427  EV_sort(zvec, wr, wi, resid, kdim);
428 
429  if (resid[nvec-1] < m_evtol)
430  {
431  idone = nvec;
432  }
433 
434  if (m_comm->GetRank() == 0)
435  {
436  evlout << "-- Iteration = " << itrn << ", H(k+1, k) = "
437  << resnorm << endl;
438  evlout.precision(4);
439  evlout.setf(ios::scientific, ios::floatfield);
441  {
442  evlout << " Magnitude Angle Growth "
443  << "Frequency Residual" << endl;
444  }
445  else
446  {
447  evlout << " Real Imaginary inverse real "
448  << "inverse imag Residual" << endl;
449  }
450 
451  for (int i = 0; i < kdim; i++)
452  {
453  WriteEvs(evlout,i,wr[i],wi[i],resid[i]);
454  }
455  }
456 
457  resid0 = resid[nvec-1];
458  return idone;
459 }
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(std::ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble, bool DumpInverse=true)
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:917
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:86

◆ EV_update()

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 305 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().

308 {
309  // Copy starting vector into first sequence element.
311  m_equ[0]->TransCoeffToPhys();
312 
313  m_equ[0]->SetTime(0.);
314  m_equ[0]->DoSolve();
315 
317  {
318  Array<OneD, MultiRegions::ExpListSharedPtr> fields;
319  fields = m_equ[0]->UpdateFields();
320 
321  //start Adjoint with latest fields of direct
322  CopyFwdToAdj();
323  m_equ[1]->TransCoeffToPhys();
324 
325  m_equ[1]->SetTime(0.);
326  m_equ[1]->DoSolve();
327  }
328 
329  // Copy starting vector into first sequence element.
331 }
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:104
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:98

◆ v_Execute()

void Nektar::SolverUtils::DriverModifiedArnoldi::v_Execute ( std::ostream &  out = std::cout)
protectedvirtual

Virtual function for solve implementation.

Implements Nektar::SolverUtils::Driver.

Reimplemented in Nektar::SolverUtils::DriverSteadyState.

Definition at line 106 of file DriverModifiedArnoldi.cpp.

References Nektar::SolverUtils::DriverArnoldi::CopyFieldToArnoldiArray(), Blas::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_imagShift, 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(), Vmath::Vcopy(), and WARNINGL0.

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

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
123  Array<OneD, NekDouble> alpha = Array<OneD, NekDouble> (m_kdim+1, 0.0);
124  Array<OneD, NekDouble> wr = Array<OneD, NekDouble> (m_kdim, 0.0);
125  Array<OneD, NekDouble> wi = Array<OneD, NekDouble> (m_kdim, 0.0);
126  Array<OneD, NekDouble> zvec = Array<OneD, NekDouble> (m_kdim*m_kdim, 0.0);
127 
128  Array<OneD, Array<OneD, NekDouble> > Kseq
129  = Array<OneD, Array<OneD, NekDouble> > (m_kdim + 1);
130  Array<OneD, Array<OneD, NekDouble> > Tseq
131  = Array<OneD, Array<OneD, NekDouble> > (m_kdim + 1);
132  for (i = 0; i < m_kdim + 1; ++i)
133  {
134  Kseq[i] = Array<OneD, NekDouble>(ntot, 0.0);
135  Tseq[i] = Array<OneD, NekDouble>(ntot, 0.0);
136  }
137 
138  // Copy starting vector into second sequence element (temporary).
139  if(m_session->DefinesFunction("InitialConditions"))
140  {
141  if (m_comm->GetRank() == 0)
142  {
143  out << "\tInital vector : specified in input file " << endl;
144  }
145  m_equ[0]->SetInitialConditions(0.0,false);
146 
147  CopyFieldToArnoldiArray(Kseq[1]);
148  }
149  else
150  {
151  if (m_comm->GetRank() == 0)
152  {
153  out << "\tInital vector : random " << endl;
154  }
155 
156  NekDouble eps=0.0001;
157  Vmath::FillWhiteNoise(ntot, eps , &Kseq[1][0], 1);
158  }
159 
160  // Perform one iteration to enforce boundary conditions.
161  // Set this as the initial value in the sequence.
162  EV_update(Kseq[1], Kseq[0]);
163  if (m_comm->GetRank() == 0)
164  {
165  out << "Iteration: " << 0 << endl;
166  }
167 
168  // Normalise first vector in sequence
169  alpha[0] = Blas::Ddot(ntot, &Kseq[0][0], 1, &Kseq[0][0], 1);
170  m_comm->AllReduce(alpha[0], Nektar::LibUtilities::ReduceSum);
171  alpha[0] = std::sqrt(alpha[0]);
172  Vmath::Smul(ntot, 1.0/alpha[0], Kseq[0], 1, Kseq[0], 1);
173 
174  // Fill initial krylov sequence
175  NekDouble resid0;
176  for (i = 1; !converged && i <= m_kdim; ++i)
177  {
178  // Compute next vector
179  EV_update(Kseq[i-1], Kseq[i]);
180 
181  // Normalise
182  alpha[i] = Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[i][0], 1);
183  m_comm->AllReduce(alpha[i], Nektar::LibUtilities::ReduceSum);
184  alpha[i] = std::sqrt(alpha[i]);
185 
186  //alpha[i] = std::sqrt(alpha[i]);
187  Vmath::Smul(ntot, 1.0/alpha[i], Kseq[i], 1, Kseq[i], 1);
188 
189  // Copy Krylov sequence into temporary storage
190  for (int k = 0; k < i + 1; ++k)
191  {
192  Vmath::Vcopy(ntot, Kseq[k], 1, Tseq[k], 1);
193  }
194 
195  // Generate Hessenberg matrix and compute eigenvalues of it.
196  EV_small(Tseq, ntot, alpha, i, zvec, wr, wi, resnorm);
197 
198  // Test for convergence.
199  converged = EV_test(i, i, zvec, wr, wi, resnorm,
200  std::min(i, m_nvec), evlout, resid0);
201  if ( i >= m_nvec)
202  {
203  converged = max (converged, 0);
204  }
205  else
206  {
207  converged = 0;
208  }
209 
210  if (m_comm->GetRank() == 0)
211  {
212  out << "Iteration: " << i << " (residual : " << resid0
213  << ")" <<endl;
214  }
215  }
216 
217  // Continue with full sequence
218  if (!converged)
219  {
220  for (i = m_kdim + 1; !converged && i <= m_nits; ++i)
221  {
222  // Shift all the vectors in the sequence.
223  // First vector is removed.
224  for (int j = 1; j <= m_kdim; ++j)
225  {
226  alpha[j-1] = alpha[j];
227  Vmath::Vcopy(ntot, Kseq[j], 1, Kseq[j-1], 1);
228  }
229 
230  // Compute next vector
231  EV_update(Kseq[m_kdim - 1], Kseq[m_kdim]);
232 
233  // Compute new scale factor
234  alpha[m_kdim] = Blas::Ddot(ntot, &Kseq[m_kdim][0], 1,
235  &Kseq[m_kdim][0], 1);
236  m_comm->AllReduce(alpha[m_kdim], Nektar::LibUtilities::ReduceSum);
237  alpha[m_kdim] = std::sqrt(alpha[m_kdim]);
238  Vmath::Smul(ntot, 1.0/alpha[m_kdim], Kseq[m_kdim], 1,
239  Kseq[m_kdim], 1);
240 
241  // Copy Krylov sequence into temporary storage
242  for (int k = 0; k < m_kdim + 1; ++k)
243  {
244  Vmath::Vcopy(ntot, Kseq[k], 1, Tseq[k], 1);
245  }
246 
247  // Generate Hessenberg matrix and compute eigenvalues of it
248  EV_small(Tseq, ntot, alpha, m_kdim, zvec, wr, wi, resnorm);
249 
250  // Test for convergence.
251  converged = EV_test(i, m_kdim, zvec, wr, wi, resnorm,
252  m_nvec, evlout, resid0);
253 
254  if (m_comm->GetRank() == 0)
255  {
256  out << "Iteration: " << i << " (residual : "
257  << resid0 << ")" <<endl;
258  }
259  }
260  }
261 
262  m_equ[0]->Output();
263 
264  // Evaluate and output computation time and solution accuracy.
265  // The specific format of the error output is essential for the
266  // regression tests to work.
267  // Evaluate L2 Error
268  for(j = 0; j < m_equ[0]->GetNvariables(); ++j)
269  {
270  NekDouble vL2Error = m_equ[0]->L2Error(j,false);
271  NekDouble vLinfError = m_equ[0]->LinfError(j);
272  if (m_comm->GetRank() == 0)
273  {
274  out << "L 2 error (variable " << m_equ[0]->GetVariable(j)
275  << ") : " << vL2Error << endl;
276  out << "L inf error (variable " << m_equ[0]->GetVariable(j)
277  << ") : " << vLinfError << endl;
278  }
279  }
280 
281  // Process eigenvectors and write out.
282  EV_post(Tseq, Kseq, ntot, min(--i, m_kdim), m_nvec, zvec, wr, wi,
283  converged);
284 
285  WARNINGL0(m_imagShift == 0,"Complex Shift applied. "
286  "Need to implement Ritz re-evaluation of"
287  "eigenvalue. Only one half of complex "
288  "value will be correct");
289 
290  // store eigenvalues so they can be accessed from driver class
291  m_real_evl = wr;
292  m_imag_evl = wi;
293 
294  // Close the runtime info file.
295  if (m_comm->GetRank() == 0)
296  {
297  evlout.close();
298  }
299 }
Array< OneD, NekDouble > m_imag_evl
Definition: DriverArnoldi.h:69
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, std::ofstream &evlout, NekDouble &resid0)
Tests for convergence of eigenvalues of H.
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:216
int m_nvec
Dimension of Krylov subspace.
Definition: DriverArnoldi.h:55
Array< OneD, NekDouble > m_real_evl
Operator in solve call is negated.
Definition: DriverArnoldi.h:68
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:223
double NekDouble
int m_nfields
interval to dump information if required.
Definition: DriverArnoldi.h:63
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:140
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:98
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:139
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:86
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:89
int m_nits
Number of vectors to test.
Definition: DriverArnoldi.h:56

◆ v_InitObject()

void Nektar::SolverUtils::DriverModifiedArnoldi::v_InitObject ( std::ostream &  out = std::cout)
protectedvirtual

Virtual function for initialisation implementation.

Reimplemented from Nektar::SolverUtils::DriverArnoldi.

Reimplemented in Nektar::SolverUtils::DriverSteadyState.

Definition at line 78 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().

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 }
SOLVER_UTILS_EXPORT void ArnoldiSummary(std::ostream &out)
virtual void v_InitObject(std::ostream &out=std::cout)
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:98
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:86
int m_nequ
number of equations
Definition: Driver.h:101

Friends And Related Function Documentation

◆ MemoryManager< DriverModifiedArnoldi >

friend class MemoryManager< DriverModifiedArnoldi >
friend

Definition at line 49 of file DriverModifiedArnoldi.h.

Member Data Documentation

◆ className

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

Name of the class.

Definition at line 63 of file DriverModifiedArnoldi.h.

◆ driverLookupId

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

Definition at line 134 of file DriverModifiedArnoldi.h.