Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
DriverArnoldi.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File DriverArnoldi.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Base Driver class for the stability solver
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 #include <iomanip>
36 
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43 namespace SolverUtils
44 {
45 
46 /**
47  * Constructor
48  */
49 DriverArnoldi::DriverArnoldi(const LibUtilities::SessionReaderSharedPtr pSession)
50  : Driver(pSession)
51 {
52  m_session->LoadParameter("IO_InfoSteps", m_infosteps, 1);
53 };
54 
55 
56 /**
57  * Destructor
58  */
60 {
61 };
62 
63 
64 /**
65  * Arnoldi driver initialisation
66  */
67 void DriverArnoldi::v_InitObject(ostream &out)
68 {
70  m_session->MatchSolverInfo("SolverType",
71  "VelocityCorrectionScheme",
73 
75  {
76  m_period = m_session->GetParameter("TimeStep")
77  * m_session->GetParameter("NumSteps");
78  m_nfields = m_equ[0]->UpdateFields().num_elements() - 1;
79 
80  }
81  else
82  {
83  m_period = 1.0;
84  m_nfields = m_equ[0]->UpdateFields().num_elements();
85  }
86 
87  if(m_session->DefinesSolverInfo("ModeType") &&
88  (boost::iequals(m_session->GetSolverInfo("ModeType"),
89  "SingleMode")||
90  boost::iequals(m_session->GetSolverInfo("ModeType"),
91  "HalfMode")))
92  {
93  for(int i = 0; i < m_nfields; ++i)
94  {
95  m_equ[0]->UpdateFields()[i]->SetWaveSpace(true);
96  }
97  }
98  m_negatedOp = m_equ[0]->v_NegatedOp();
99 
100  m_session->LoadParameter("kdim", m_kdim, 16);
101  m_session->LoadParameter("nvec", m_nvec, 2);
102  m_session->LoadParameter("nits", m_nits, 500);
103  m_session->LoadParameter("evtol", m_evtol, 1e-06);
104 
105  ASSERTL0( m_kdim >= m_nvec, "nvec cannot be larger than kdim.");
106 
107  m_session->LoadParameter("realShift", m_realShift, 0.0);
108  m_equ[0]->SetLambda(m_realShift);
109 
110  m_session->LoadParameter("imagShift", m_imagShift, 0.0);
111 
112 }
113 
114 void DriverArnoldi::ArnoldiSummary(std::ostream &out)
115 {
116  if (m_comm->GetRank() == 0)
117  {
118  if(m_session->DefinesSolverInfo("ModeType") &&
119  boost::iequals(m_session->GetSolverInfo("ModeType"),
120  "SingleMode"))
121  {
122  out << "\tSingle Fourier mode : true " << endl;
123  ASSERTL0(m_session->DefinesSolverInfo("Homogeneous"),
124  "Expected a homogeneous expansion to be defined "
125  "with single mode");
126  }
127  else
128  {
129  out << "\tSingle Fourier mode : false " << endl;
130  }
131  if(m_session->DefinesSolverInfo("BetaZero"))
132  {
133  out << "\tBeta set to Zero : true (overrides LHom)"
134  << endl;
135  }
136  else
137  {
138  out << "\tBeta set to Zero : false " << endl;
139  }
140 
142  {
143  out << "\tEvolution operator : "
144  << m_session->GetSolverInfo("EvolutionOperator")
145  << endl;
146  }
147  else
148  {
149  out << "\tShift (Real,Imag) : " << m_realShift
150  << "," << m_imagShift << endl;
151  }
152  out << "\tKrylov-space dimension : " << m_kdim << endl;
153  out << "\tNumber of vectors : " << m_nvec << endl;
154  out << "\tMax iterations : " << m_nits << endl;
155  out << "\tEigenvalue tolerance : " << m_evtol << endl;
156  out << "======================================================"
157  << endl;
158  }
159 }
160 
161 /**
162  * Copy Arnoldi array to field variables which depend from
163  * either the m_fields or m_forces
164  */
166 {
167 
168  Array<OneD, MultiRegions::ExpListSharedPtr>& fields = m_equ[0]->UpdateFields();
169  int nq = fields[0]->GetNcoeffs();
170 
171  for (int k = 0; k < m_nfields; ++k)
172  {
173  Vmath::Vcopy(nq, &array[k*nq], 1, &fields[k]->UpdateCoeffs()[0], 1);
174  fields[k]->SetPhysState(false);
175  }
176 };
177 
178 /**
179  * Copy field variables which depend from either the m_fields
180  * or m_forces array the Arnoldi array
181  */
183 {
184 
186 
188  {
189  // This matters for the Adaptive SFD method because
190  // m_equ[1] is the nonlinear problem with non
191  // homogeneous BCs.
192  fields = m_equ[0]->UpdateFields();
193  }
194  else
195  {
196  fields = m_equ[m_nequ-1]->UpdateFields();
197  }
198 
199  for (int k = 0; k < m_nfields; ++k)
200  {
201  int nq = fields[0]->GetNcoeffs();
202  Vmath::Vcopy(nq, &fields[k]->GetCoeffs()[0], 1, &array[k*nq], 1);
203  fields[k]->SetPhysState(false);
204 
205  }
206 };
207 
208 
209 /**
210  * Initialisation for the transient growth
211  */
213 {
215 
217  {
218  fields = m_equ[0]->UpdateFields();
219  int nq = fields[0]->GetNcoeffs();
220 
221  for (int k=0 ; k < m_nfields; ++k)
222  {
223  Vmath::Vcopy(nq,
224  &fields[k]->GetCoeffs()[0], 1,
225  &m_equ[1]->UpdateFields()[k]->UpdateCoeffs()[0], 1);
226 
227  }
228  }
229  else
230  {
231  ASSERTL0(false,"Transient Growth non available for Coupled Solver");
232 
233  }
234 };
235 
236 void DriverArnoldi::WriteFld(std::string file, std::vector<Array<OneD, NekDouble> > coeffs)
237 {
238 
239  std::vector<std::string> variables(m_nfields);
240 
241  ASSERTL1(coeffs.size() >= m_nfields, "coeffs is not of the correct length");
242  for(int i = 0; i < m_nfields; ++i)
243  {
244  variables[i] = m_equ[0]->GetVariable(i);
245  }
246 
247  m_equ[0]->WriteFld(file,m_equ[0]->UpdateFields()[0], coeffs, variables);
248 }
249 
250 
251 void DriverArnoldi::WriteFld(std::string file, Array<OneD, NekDouble> coeffs)
252 {
253 
254  std::vector<std::string> variables(m_nfields);
255  std::vector<Array<OneD, NekDouble> > fieldcoeffs(m_nfields);
256 
257  int ncoeffs = m_equ[0]->UpdateFields()[0]->GetNcoeffs();
258  ASSERTL1(coeffs.num_elements() >= ncoeffs*m_nfields,"coeffs is not of sufficient size");
259 
260  for(int i = 0; i < m_nfields; ++i)
261  {
262  variables[i] = m_equ[0]->GetVariable(i);
263  fieldcoeffs[i] = coeffs + i*ncoeffs;
264  }
265 
266  m_equ[0]->WriteFld(file,m_equ[0]->UpdateFields()[0], fieldcoeffs, variables);
267 }
268 
270  ostream &evlout,
271  const int i,
272  const NekDouble re_ev,
273  const NekDouble im_ev,
274  NekDouble resid,
275  bool DumpInverse)
276 {
278  {
279  NekDouble abs_ev = hypot (re_ev, im_ev);
280  NekDouble ang_ev = atan2 (im_ev, re_ev);
281 
282  evlout << "EV: " << setw(2) << i
283  << setw(12) << abs_ev
284  << setw(12) << ang_ev
285  << setw(12) << log (abs_ev) / m_period
286  << setw(12) << ang_ev / m_period;
287 
288  if(resid != NekConstants::kNekUnsetDouble)
289  {
290  evlout << setw(12) << resid;
291  }
292  evlout << endl;
293  }
294  else
295  {
296  NekDouble invmag = 1.0/(re_ev*re_ev + im_ev*im_ev);
297  NekDouble sign;
298  if(m_negatedOp)
299  {
300  sign = -1.0;
301  }
302  else
303  {
304  sign = 1.0;
305  }
306 
307  evlout << "EV: " << setw(2) << i
308  << setw(14) << sign*re_ev
309  << setw(14) << sign*im_ev;
310 
311  if(DumpInverse)
312  {
313  evlout << setw(14) << sign*re_ev*invmag + m_realShift
314  << setw(14) << sign*im_ev*invmag + m_imagShift;
315  }
316 
317  if(resid != NekConstants::kNekUnsetDouble)
318  {
319  evlout << setw(12) << resid;
320  }
321  evlout << endl;
322  }
323 }
324 
325 }
326 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
void CopyFwdToAdj()
Copy the forward field to the adjoint system in transient growth calculations.
NekDouble m_evtol
Maxmum number of iterations.
Definition: DriverArnoldi.h:58
void CopyArnoldiArrayToField(Array< OneD, NekDouble > &array)
Copy Arnoldi storage to fields.
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:22
SOLVER_UTILS_EXPORT void ArnoldiSummary(std::ostream &out)
STL namespace.
void CopyFieldToArnoldiArray(Array< OneD, NekDouble > &array)
Copy fields to Arnoldi storage.
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
void WriteFld(std::string file, std::vector< Array< OneD, NekDouble > > coeffs)
Write coefficients to file.
virtual SOLVER_UTILS_EXPORT void v_InitObject(std::ostream &out=std::cout)
Definition: Driver.cpp:92
virtual ~DriverArnoldi()
Destructor.
int m_nvec
Dimension of Krylov subspace.
Definition: DriverArnoldi.h:56
void WriteEvs(std::ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble, bool DumpInverse=true)
virtual void v_InitObject(std::ostream &out=std::cout)
enum EvolutionOperatorType m_EvolutionOperator
Evolution Operator.
Definition: Driver.h:100
double NekDouble
static const NekDouble kNekUnsetDouble
int m_nfields
interval to dump information if required.
Definition: DriverArnoldi.h:64
bool m_timeSteppingAlgorithm
Period of time stepping algorithm.
Definition: DriverArnoldi.h:60
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:94
NekDouble m_period
Tolerance of iteratiosn.
Definition: DriverArnoldi.h:59
int m_infosteps
underlying operator is time stepping
Definition: DriverArnoldi.h:62
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:85
Base class for the development of solvers.
Definition: Driver.h:65
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
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_nequ
number of equations
Definition: Driver.h:97
int m_nits
Number of vectors to test.
Definition: DriverArnoldi.h:57