Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 (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)
 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)
 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,
EquationSystemSharedPtr
GetEqu ()
 
SOLVER_UTILS_EXPORT Array
< OneD, NekDouble
GetRealEvl (void)
 
SOLVER_UTILS_EXPORT Array
< OneD, NekDouble
GetImagEvl (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...
 
Array< OneD,
EquationSystemSharedPtr
m_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 58 of file DriverModifiedArnoldi.cpp.

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

Destructor.

Definition at line 68 of file DriverModifiedArnoldi.cpp.

69 {
70 }

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(), and CellMLToNektar.cellml_metadata::p.

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

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

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

557 {
558  NekDouble wgt, norm;
559 
560  // Generate big eigenvectors
561  for (int j = 0; j < nvec; ++j)
562  {
563  Vmath::Zero(ntot, evecs[j], 1);
564  for (int i = 0; i < kdim; ++i)
565  {
566  wgt = zvec[i + j*kdim];
567  Vmath::Svtvp(ntot, wgt, bvecs[i], 1, evecs[j], 1, evecs[j], 1);
568  }
569  }
570 
571  // Normalise the big eigenvectors
572  for (int i = 0; i < nvec; ++i)
573  {
574  if (wi[i] == 0.0) // Real mode
575  {
576  norm = Blas::Ddot(ntot, &evecs[i][0], 1, &evecs[i][0], 1);
577  m_comm->AllReduce(norm, Nektar::LibUtilities::ReduceSum);
578  norm = std::sqrt(norm);
579  Vmath::Smul(ntot, 1.0/norm, evecs[i], 1, evecs[i], 1);
580  }
581  else
582  {
583  norm = Blas::Ddot(ntot, &evecs[i][0], 1, &evecs[i][0], 1);
584  norm += Blas::Ddot(ntot, &evecs[i+1][0], 1, &evecs[i+1][0], 1);
585  m_comm->AllReduce(norm, Nektar::LibUtilities::ReduceSum);
586  norm = std::sqrt(norm);
587 
588  Vmath::Smul(ntot, 1.0/norm, evecs[i], 1, evecs[i], 1);
589  Vmath::Smul(ntot, 1.0/norm, evecs[i+1], 1, evecs[i+1], 1);
590 
591  i++;
592  }
593  }
594 }
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:485
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:213
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:436
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:85
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:373
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 498 of file DriverModifiedArnoldi.cpp.

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

Referenced by v_Execute().

508 {
509  if (icon == 0)
510  {
511  // Not converged, write final Krylov vector
512  ASSERTL0(false, "Convergence was not achieved within the "
513  "prescribed number of iterations.");
514  }
515  else if (icon < 0)
516  {
517  // Minimum residual reached
518  ASSERTL0(false, "Minimum residual reached.");
519  }
520  else if (icon == nvec)
521  {
522  // Converged, write out eigenvectors
523  EV_big(Tseq, Kseq, ntot, kdim, icon, zvec, wr, wi);
524  Array<OneD, MultiRegions::ExpListSharedPtr> fields
525  = m_equ[0]->UpdateFields();
526 
527  for (int j = 0; j < icon; ++j)
528  {
529  std::string file = m_session->GetSessionName() + "_eig_"
530  + boost::lexical_cast<std::string>(j)
531  + ".fld";
532 
533  WriteEvs(cout, j, wr[j], wi[j]);
534  WriteFld(file,Kseq[j]);
535  }
536  }
537  else
538  {
539  // Not recognised
540  ASSERTL0(false, "Unrecognised value.");
541  }
542 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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: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.

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

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

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

Sorts the computed eigenvalues by smallest residual first.

Definition at line 461 of file DriverModifiedArnoldi.cpp.

References Vmath::Vcopy().

Referenced by EV_test().

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

407 {
408  int idone = 0;
409  // NekDouble period = 0.1;
410 
411  Array<OneD, NekDouble> resid(kdim);
412  for (int i = 0; i < kdim; ++i)
413  {
414  NekDouble tmp = std::sqrt(Vmath::Dot(kdim, &zvec[0] + i*kdim, 1,
415  &zvec[0] + i*kdim, 1));
416  resid[i] = resnorm * std::fabs(zvec[kdim - 1 + i*kdim]) / tmp;
417  if (wi[i] < 0.0)
418  {
419  resid[i-1] = resid[i] = hypot(resid[i-1], resid[i]);
420  }
421  }
422 
423  EV_sort(zvec, wr, wi, resid, kdim);
424 
425  if (resid[nvec-1] < m_evtol)
426  {
427  idone = nvec;
428  }
429 
430  if (m_comm->GetRank() == 0)
431  {
432  evlout << "-- Iteration = " << itrn << ", H(k+1, k) = "
433  << resnorm << endl;
434  evlout.precision(4);
435  evlout.setf(ios::scientific, ios::floatfield);
437  {
438  evlout << " Magnitude Angle Growth "
439  << "Frequency Residual" << endl;
440  }
441  else
442  {
443  evlout << " Real Imaginary inverse real "
444  << "inverse imag Residual" << endl;
445  }
446 
447  for (int i = 0; i < kdim; i++)
448  {
449  WriteEvs(evlout,i,wr[i],wi[i],resid[i]);
450  }
451  }
452 
453  resid0 = resid[nvec-1];
454  return idone;
455 }
NekDouble m_evtol
Maxmum number of iterations.
Definition: DriverArnoldi.h:58
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:60
T Dot(int n, const T *w, const T *x)
vvtvp (vector times vector times vector): z = w*x*y
Definition: Vmath.cpp:914
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 303 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().

306 {
307  // Copy starting vector into first sequence element.
309  m_equ[0]->TransCoeffToPhys();
310 
311  m_equ[0]->DoSolve();
312 
314  {
315  Array<OneD, MultiRegions::ExpListSharedPtr> fields;
316  fields = m_equ[0]->UpdateFields();
317 
318  //start Adjoint with latest fields of direct
319  CopyFwdToAdj();
320  m_equ[1]->TransCoeffToPhys();
321 
322  m_equ[1]->DoSolve();
323  }
324 
325  // Copy starting vector into first sequence element.
327 }
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 ( std::ostream &  out = std::cout)
protectedvirtual

Virtual function for solve implementation.

Implements Nektar::SolverUtils::Driver.

Reimplemented in Nektar::SolverUtils::DriverSteadyState.

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

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

77 {
79 
80  m_equ[0]->PrintSummary(out);
81 
82  // Print session parameters
83  if (m_comm->GetRank() == 0)
84  {
85  out << "\tArnoldi solver type : Modified Arnoldi" << endl;
86  }
87 
89 
90  for( int i = 0; i < m_nequ; ++i)
91  {
92  m_equ[i]->DoInitialise();
93  }
94 
95  //FwdTrans Initial conditions to be in Coefficient Space
96  m_equ[m_nequ-1] ->TransPhysToCoeff();
97 
98 }
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:94
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:85
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.