Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator 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().

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

554 {
555  NekDouble wgt, norm;
556 
557  // Generate big eigenvectors
558  for (int j = 0; j < nvec; ++j)
559  {
560  Vmath::Zero(ntot, evecs[j], 1);
561  for (int i = 0; i < kdim; ++i)
562  {
563  wgt = zvec[i + j*kdim];
564  Vmath::Svtvp(ntot, wgt, bvecs[i], 1, evecs[j], 1, evecs[j], 1);
565  }
566  }
567 
568  // Normalise the big eigenvectors
569  for (int i = 0; i < nvec; ++i)
570  {
571  if (wi[i] == 0.0) // Real mode
572  {
573  norm = Blas::Ddot(ntot, &evecs[i][0], 1, &evecs[i][0], 1);
574  m_comm->AllReduce(norm, Nektar::LibUtilities::ReduceSum);
575  norm = std::sqrt(norm);
576  Vmath::Smul(ntot, 1.0/norm, evecs[i], 1, evecs[i], 1);
577  }
578  else
579  {
580  norm = Blas::Ddot(ntot, &evecs[i][0], 1, &evecs[i][0], 1);
581  norm += Blas::Ddot(ntot, &evecs[i+1][0], 1, &evecs[i+1][0], 1);
582  m_comm->AllReduce(norm, Nektar::LibUtilities::ReduceSum);
583  norm = std::sqrt(norm);
584 
585  Vmath::Smul(ntot, 1.0/norm, evecs[i], 1, evecs[i], 1);
586  Vmath::Smul(ntot, 1.0/norm, evecs[i+1], 1, evecs[i+1], 1);
587 
588  i++;
589  }
590  }
591 }
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:434
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

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

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

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

References Vmath::Vcopy().

Referenced by EV_test().

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

404 {
405  int idone = 0;
406  // NekDouble period = 0.1;
407 
408  Array<OneD, NekDouble> resid(kdim);
409  for (int i = 0; i < kdim; ++i)
410  {
411  NekDouble tmp = std::sqrt(Vmath::Dot(kdim, &zvec[0] + i*kdim, 1,
412  &zvec[0] + i*kdim, 1));
413  resid[i] = resnorm * std::fabs(zvec[kdim - 1 + i*kdim]) / tmp;
414  if (wi[i] < 0.0)
415  {
416  resid[i-1] = resid[i] = hypot(resid[i-1], resid[i]);
417  }
418  }
419 
420  EV_sort(zvec, wr, wi, resid, kdim);
421 
422  if (resid[nvec-1] < m_evtol)
423  {
424  idone = nvec;
425  }
426 
427  if (m_comm->GetRank() == 0)
428  {
429  evlout << "-- Iteration = " << itrn << ", H(k+1, k) = "
430  << resnorm << endl;
431  evlout.precision(4);
432  evlout.setf(ios::scientific, ios::floatfield);
434  {
435  evlout << " Magnitude Angle Growth "
436  << "Frequency Residual" << endl;
437  }
438  else
439  {
440  evlout << " Real Imaginary inverse real "
441  << "inverse imag Residual" << endl;
442  }
443 
444  for (int i = 0; i < kdim; i++)
445  {
446  WriteEvs(evlout,i,wr[i],wi[i],resid[i]);
447  }
448  }
449 
450  resid0 = resid[nvec-1];
451  return idone;
452 }
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:900
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 300 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().

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

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