Nektar++
DriverArpack.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: DriverArpack.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// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Arnoldi solver using Arpack
32//
33///////////////////////////////////////////////////////////////////////////////
34
36
37using namespace std;
38
39namespace Nektar::SolverUtils
40{
41
44 "LargestReal", 0),
46 "SmallestReal", 1),
48 "LargestImag", 2),
50 "SmallestImag", 3),
52 "LargestMag", 4),
54 "SmallestMag", 5),
55};
58 "LargestMag");
61
62std::string DriverArpack::className =
64
65std::string DriverArpack::ArpackProblemTypeTrans[6] = {"LR", "SR", "LI",
66 "SI", "LM", "SM"};
67
68/**
69 *
70 */
73 : DriverArnoldi(pSession, pGraph)
74{
75}
76
77/**
78 *
79 */
81{
83
84 // Initialisation of Arnoldi parameters
85 m_maxn = 1000000; // Maximum size of the problem
86 m_maxnev = 200; // maximum number of eigenvalues requested
87 m_maxncv = 500; // Largest number of basis vector used in Implicitly
88 // Restarted Arnoldi
89
90 // Error alerts
91 ASSERTL0(m_nvec <= m_maxnev, "NEV is greater than MAXNEV");
92 ASSERTL0(m_kdim <= m_maxncv, "NEV is greater than MAXNEV");
93 ASSERTL0(2 <= m_kdim - m_nvec, "NCV-NEV is less than 2");
94
95 ASSERTL0(m_comm->GetSize() == 1,
96 "..ARPACK is not currently set-up for parallel execution...\n");
97
98 m_equ[0]->PrintSummary(out);
99
100 // Print session parameters
101 out << "\tArnoldi solver type : Arpack" << endl;
102 out << "\tArpack problem type : ";
103 out << ArpackProblemTypeTrans[m_session->GetSolverInfoAsEnum<int>(
104 "ArpackProblemType")]
105 << endl;
107
108 // Initialization
109 for (int i = 0; i < m_nequ; ++i)
110 {
111 m_equ[i]->DoInitialise();
112 }
113
114 // FwdTrans Initial conditions to be in Coefficient Space
115 m_equ[m_nequ - 1]->TransPhysToCoeff();
116}
117
118/**
119 *
120 */
121void DriverArpack::v_Execute(ostream &out)
122
123{
124 Array<OneD, NekDouble> tmpworkd;
125
126 int nq = m_equ[0]
127 ->UpdateFields()[0]
128 ->GetNcoeffs(); // Number of points in the mesh
129 int n = m_nfields * nq; // Number of points in eigenvalue calculation
130 int lworkl = 3 * m_kdim * (m_kdim + 2); // Size of work array
131 int ido; // REVERSE COMMUNICATION parameter. At the first call must be
132 // initialised at 0
133 int info; // do not set initial vector (info=0 random initial vector, info=1
134 // read initial vector from session file)
135
136 int iparam[11];
137 int ipntr[14];
138
139 Array<OneD, int> ritzSelect;
144 NekDouble sigmar, sigmai;
145
146 Array<OneD, NekDouble> resid(n);
148 Array<OneD, NekDouble> workl(lworkl, 0.0);
149 Array<OneD, NekDouble> workd(3 * n, 0.0);
150
151 ASSERTL0(n <= m_maxn, "N is greater than MAXN");
152
153 if (m_session->DefinesFunction("InitialConditions"))
154 {
155 out << "\tInital vector : input file " << endl;
156 info = 1;
158 }
159 else
160 {
161 out << "\tInital vector : random " << endl;
162 info = 0;
163 }
164
165 char B;
166
167 iparam[0] = 1; // strategy for shift-invert
168 iparam[1] = 0; // (deprecated)
169 iparam[2] = m_nits; // maximum number of iterations allowed/taken
170 iparam[3] = 1; // blocksize to be used for recurrence
171 iparam[4] = 0; // number of converged ritz eigenvalues
172 iparam[5] = 0; // (deprecated)
173
174 // Use generalized B matrix for coupled solver.
176 {
177 iparam[6] = 1; // computation mode 1=> matrix-vector prod
178 B = 'I';
179 }
180 else
181 {
182 iparam[6] = 3; // computation mode 1=> matrix-vector prod
183 B = 'G';
184 }
185#if 0
186 if((fabs(m_realShift) > NekConstants::kNekZeroTol)|| // use shift if m_realShift > 1e-12
188 {
189 iparam[6] = 3; // This was 3 need to know what to set it to
190 B = 'G';
191 }
192 else
193 {
194 iparam[6] = 1; // computation mode 1=> matrix-vector prod
195 B = 'I';
196 }
197#endif
198 iparam[7] = 0; // (for shift-invert)
199 iparam[8] = 0; // number of MV operations
200 iparam[9] = 0; // number of BV operations
201 iparam[10] = 0; // number of reorthogonalisation steps
202
203 int cycle = 0;
204 const char *problem =
205 ArpackProblemTypeTrans[m_session->GetSolverInfoAsEnum<int>(
206 "ArpackProblemType")]
207 .c_str();
208
209 std::string name = m_session->GetSessionName() + ".evl";
210 ofstream pFile(name.c_str());
211
212 ido = 0; // At the first call must be initialisedat 0
213
214 while (ido != 99) // ido==-1 || ido==1 || ido==0)
215 {
216 // Routine for eigenvalue evaluation for non-symmetric operators
217 Arpack::Dnaupd(ido, &B, // B='I' for std eval problem
218 n, problem, m_nvec, m_evtol, &resid[0], m_kdim, &v[0], n,
219 iparam, ipntr, &workd[0], &workl[0], lworkl, info);
220
221 // Plotting of real and imaginary part of the
222 // eigenvalues from workl
223 out << "\rIteration " << cycle << ", output: " << info
224 << ", ido=" << ido << " " << std::flush;
225
226 if (!((cycle - 1) % m_kdim) && (cycle > m_kdim) && (ido != 2))
227 {
228 pFile << "Krylov spectrum at iteration: " << cycle << endl;
229
231 {
232 pFile << "EV Magnitude Angle Growth Frequency "
233 "Residual"
234 << endl;
235 }
236 else
237 {
238 pFile << "EV Real Imaginary inverse real inverse "
239 "imag Residual"
240 << endl;
241 }
242
243 out << endl;
244 for (int k = 0; k < m_kdim; ++k)
245 {
246 // write m_kdim eigs to screen
247 WriteEvs(pFile, k, workl[ipntr[5] - 1 + k],
248 workl[ipntr[6] - 1 + k]);
249 }
250 }
251
252 if (ido == 99)
253 {
254 break;
255 }
256
257 switch (ido)
258 {
259 case -1:
260 case 1: // Note that ido=1 we are using input x
261 // (workd[inptr[0]-1]) rather than Mx as
262 // recommended in manual since it is not
263 // possible to impose forcing directly.
264 CopyArnoldiArrayToField(tmpworkd = workd + (ipntr[0] - 1));
265
266 m_equ[0]->TransCoeffToPhys();
267
268 m_equ[0]->DoSolve();
270 {
271 // start Adjoint with latest fields of direct
272 CopyFwdToAdj();
273
274 m_equ[1]->TransCoeffToPhys();
275 m_equ[1]->DoSolve();
276 }
277
278 if (!(cycle % m_infosteps))
279 {
280 out << endl;
281 m_equ[0]->Output();
282 }
283
284 // operated fields are copied into workd[inptr[1]-1]
285 CopyFieldToArnoldiArray(tmpworkd = workd + (ipntr[1] - 1));
286
287 cycle++;
288 break;
289 case 2: // provide y = M x (bwd trans and iproduct);
290 {
291 // workd[inptr[0]-1] copied into operator fields
292 CopyArnoldiArrayToField(tmpworkd = workd + (ipntr[0] - 1));
293
294 m_equ[0]->TransCoeffToPhys();
295
297 m_equ[0]->UpdateFields();
298 for (int i = 0; i < fields.size(); ++i)
299 {
300 fields[i]->IProductWRTBase(fields[i]->GetPhys(),
301 fields[i]->UpdateCoeffs());
302 }
303
304 // operated fields are copied into workd[inptr[1]-1]
305 CopyFieldToArnoldiArray(tmpworkd = workd + (ipntr[1] - 1));
306 break;
307 }
308 default:
309 ASSERTL0(false, "Unexpected reverse communication request.");
310 }
311 }
312
313 out << endl << "Converged in " << iparam[8] << " iterations" << endl;
314
315 ASSERTL0(info >= 0, " Error with Dnaupd");
316
317 ritzSelect = Array<OneD, int>(m_kdim, 0);
318 dr = Array<OneD, NekDouble>(m_nvec + 1, 0.0);
319 di = Array<OneD, NekDouble>(m_nvec + 1, 0.0);
320 workev = Array<OneD, NekDouble>(3 * m_kdim);
321 z = Array<OneD, NekDouble>(n * (m_nvec + 1));
322
323 if (m_negatedOp)
324 {
325 sigmar = -m_realShift;
326 }
327 else
328 {
329 sigmar = m_realShift;
330 }
331
332 // Do not pass imaginary shift to Arpack since we have not
333 // used a Fortran complex number format and so processing
334 // is mucked up. Need to do some processing afterwards.
335 sigmai = 0;
336
337 // Setting 'A', Ritz vectors are computed. 'S' for Shur vectors
338 Arpack::Dneupd(1, "A", ritzSelect.get(), dr.get(), di.get(), z.get(), n,
339 sigmar, sigmai, workev.get(), &B, n, problem, m_nvec,
340 m_evtol, resid.get(), m_kdim, v.get(), n, iparam, ipntr,
341 workd.get(), workl.get(), lworkl, info);
342
343 ASSERTL0(info == 0, " Error with Dneupd");
344
345 int nconv = iparam[4];
346
347 // Subtract off complex shift if it exists
348 if (m_negatedOp)
349 {
350 Vmath::Sadd(nconv, m_imagShift, di, 1, di, 1);
351 }
352 else
353 {
354 Vmath::Sadd(nconv, -m_imagShift, di, 1, di, 1);
355 }
356
357 WARNINGL0(m_imagShift == 0, "Complex Shift applied. "
358 "Need to implement Ritz re-evaluation of"
359 "eigenvalue. Only one half of complex "
360 "value will be correct");
361
363 m_equ[0]->UpdateFields();
364
365 out << "Converged Eigenvalues: " << nconv << endl;
366 pFile << "Converged Eigenvalues: " << nconv << endl;
367
369 {
370 out << " Magnitude Angle Growth Frequency" << endl;
371 pFile << " Magnitude Angle Growth Frequency"
372 << endl;
373 for (int i = 0; i < nconv; ++i)
374 {
375 WriteEvs(out, i, dr[i], di[i]);
376 WriteEvs(pFile, i, dr[i], di[i]);
377
378 std::string file = m_session->GetSessionName() + "_eig_" +
379 std::to_string(i) + ".fld";
380 WriteFld(file, z + i * n);
381 }
382 }
383 else
384 {
385 out << " Real Imaginary " << endl;
386 pFile << " Real Imaginary " << endl;
387 for (int i = 0; i < nconv; ++i)
388 {
389 WriteEvs(out, i, dr[i], di[i], NekConstants::kNekUnsetDouble,
390 false);
391 WriteEvs(pFile, i, dr[i], di[i], NekConstants::kNekUnsetDouble,
392 false);
393
394 std::string file = m_session->GetSessionName() + "_eig_" +
395 std::to_string(i) + ".fld";
396 WriteFld(file, z + i * n);
397 }
398 }
399
400 m_real_evl = dr;
401 m_imag_evl = di;
402
403 pFile.close();
404
405 for (int j = 0; j < m_equ[0]->GetNvariables(); ++j)
406 {
407 NekDouble vL2Error = m_equ[0]->L2Error(j, false);
408 NekDouble vLinfError = m_equ[0]->LinfError(j);
409 if (m_comm->GetRank() == 0)
410 {
411 out << "L 2 error (variable " << m_equ[0]->GetVariable(j)
412 << ") : " << vL2Error << endl;
413 out << "L inf error (variable " << m_equ[0]->GetVariable(j)
414 << ") : " << vLinfError << endl;
415 }
416 }
417}
418
419} // namespace Nektar::SolverUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:215
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
static std::string RegisterDefaultSolverInfo(const std::string &pName, const std::string &pValue)
Registers the default string value of a solver info property.
Base class for the development of solvers.
Definition: DriverArnoldi.h:45
void CopyFwdToAdj()
Copy the forward field to the adjoint system in transient growth calculations.
void v_InitObject(std::ostream &out=std::cout) override
Virtual function for initialisation implementation.
void WriteFld(std::string file, std::vector< Array< OneD, NekDouble > > coeffs)
Write coefficients to file.
int m_infosteps
underlying operator is time stepping
Definition: DriverArnoldi.h:67
void CopyFieldToArnoldiArray(Array< OneD, NekDouble > &array)
Copy fields to Arnoldi storage.
int m_nvec
Dimension of Krylov subspace.
Definition: DriverArnoldi.h:61
bool m_timeSteppingAlgorithm
Period of time stepping algorithm.
Definition: DriverArnoldi.h:65
int m_nits
Number of vectors to test.
Definition: DriverArnoldi.h:62
Array< OneD, NekDouble > m_imag_evl
Definition: DriverArnoldi.h:75
void CopyArnoldiArrayToField(Array< OneD, NekDouble > &array)
Copy Arnoldi storage to fields.
NekDouble m_evtol
Maxmum number of iterations.
Definition: DriverArnoldi.h:63
int m_nfields
interval to dump information if required.
Definition: DriverArnoldi.h:69
SOLVER_UTILS_EXPORT void ArnoldiSummary(std::ostream &out)
Array< OneD, NekDouble > m_real_evl
Operator in solve call is negated.
Definition: DriverArnoldi.h:74
void WriteEvs(std::ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble, bool DumpInverse=true)
static std::string arpackProblemTypeDefault
Definition: DriverArpack.h:88
static DriverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Definition: DriverArpack.h:52
DriverArpack(const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
Constructor.
void v_InitObject(std::ostream &out=std::cout) override
Virtual function for initialisation implementation.
static std::string driverLookupId
Definition: DriverArpack.h:84
static std::string arpackProblemTypeLookupIds[]
Definition: DriverArpack.h:87
void v_Execute(std::ostream &out=std::cout) override
Virtual function for solve implementation.
static std::string className
Name of the class.
Definition: DriverArpack.h:63
static std::string ArpackProblemTypeTrans[]
Definition: DriverArpack.h:89
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:83
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:80
enum EvolutionOperatorType m_EvolutionOperator
Evolution Operator.
Definition: Driver.h:98
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:92
int m_nequ
number of equations
Definition: Driver.h:95
static void Dnaupd(int &ido, const char *bmat, const int &n, const char *which, const int &nev, const double &tol, double *resid, const int &ncv, double *v, const int &ldv, int *iparam, int *ipntr, double *workd, double *workl, const int &lworkl, int &info)
Top level reverse communication interface to solve real double-precision non-symmetric problems.
Definition: Arpack.hpp:114
static void Dneupd(const int &rvec, const char *howmny, const int *select, double *dr, double *di, double *z, const int &ldz, const double &sigmar, const double &sigmai, double *workev, const char *bmat, const int &n, const char *which, const int &nev, const double &tol, double *resid, const int &ncv, double *v, const int &ldv, int *iparam, int *ipntr, double *workd, double *workl, const int &lworkl, int &info)
Post-processing routine to computed eigenvector of computed eigenvalues in Dnaupd.
Definition: Arpack.hpp:128
std::shared_ptr< SessionReader > SessionReaderSharedPtr
static const NekDouble kNekUnsetDouble
static const NekDouble kNekZeroTol
DriverFactory & GetDriverFactory()
Definition: Driver.cpp:65
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::vector< double > z(NPUPPER)
double NekDouble
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.hpp:194