Bridge++  Ver. 2.0.2
fft_xyz_1dim.cpp
Go to the documentation of this file.
1 
14 #ifdef USE_FFTWLIB
15 
16 #include "fft_xyz_1dim.h"
17 
18 #ifdef USE_FACTORY_AUTOREGISTER
19 namespace {
20  bool init = FFT_xyz_1dim::register_factory();
21 }
22 #endif
23 
24 const std::string FFT_xyz_1dim::class_name = "FFT_xyz_1dim";
25 
26 //====================================================================
27 void FFT_xyz_1dim::set_parameters(const Parameters& params)
28 {
29  std::string vlevel;
30  if (!params.fetch_string("verbose_level", vlevel)) {
31  m_vl = vout.set_verbose_level(vlevel);
32  }
33 
34  //- fetch and check input parameters
35  std::string str_fft_direction;
36 
37  int err = 0;
38  err += params.fetch_string("FFT_direction", str_fft_direction);
39 
40  if (err) {
41  vout.crucial(m_vl, "Error at %s: input parameter not found.\n", class_name.c_str());
42  exit(EXIT_FAILURE);
43  }
44 
45  set_parameters(str_fft_direction);
46 }
47 
48 
49 //====================================================================
50 void FFT_xyz_1dim::get_parameters(Parameters& params) const
51 {
52  params.set_string("FFT_direction", m_is_forward ? "Forward" : "Backward");
53 
54  params.set_string("verbose_level", vout.get_verbose_level(m_vl));
55 }
56 
57 
58 //====================================================================
59 void FFT_xyz_1dim::set_parameters(const std::string& str_fft_direction)
60 {
61  //- print input parameters
62  vout.general(m_vl, "%s:\n", class_name.c_str());
63  vout.general(m_vl, " FFT_direction = %s\n", str_fft_direction.c_str());
64 
65  //- range check
66 
67  //- store values
68  if (str_fft_direction == "Forward") {
69  m_is_forward = true;
70  } else if (str_fft_direction == "Backward") {
71  m_is_forward = false;
72  } else {
73  vout.crucial(m_vl, "Error at %s: unsupported FFT direction \"%s\"\n", class_name.c_str(), str_fft_direction.c_str());
74  exit(EXIT_FAILURE);
75  }
76 }
77 
78 
79 //====================================================================
80 void FFT_xyz_1dim::fft(Field& field)
81 {
82  const int Ndim = CommonParameters::Ndim();
83 
84  //- global lattice size
85  std::vector<int> Lsize(Ndim - 1);
86 
87  Lsize[0] = CommonParameters::Lx();
88  Lsize[1] = CommonParameters::Ly();
89  Lsize[2] = CommonParameters::Lz();
90 
91  const int Lt = CommonParameters::Lt();
92  const int Lxyz = Lsize[0] * Lsize[1] * Lsize[2];
93 
94  //- local lattice size
95  std::vector<int> Nsize(Ndim - 1);
96  Nsize[0] = CommonParameters::Nx();
97  Nsize[1] = CommonParameters::Ny();
98  Nsize[2] = CommonParameters::Nz();
99 
100  const int Nin = field.nin();
101  const int Nex = field.nex();
102 
103  std::vector<int> NPE(Ndim);
104  NPE[0] = CommonParameters::NPEx();
105  NPE[1] = CommonParameters::NPEy();
106  NPE[2] = CommonParameters::NPEz();
107  NPE[3] = CommonParameters::NPEt();
108 
109  const int NPE_xyt = NPE[0] * NPE[1] * NPE[3];
110  const int NPE_xzt = NPE[0] * NPE[2] * NPE[3];
111  const int NPE_yzt = NPE[1] * NPE[2] * NPE[3];
112 
113  if ((NPE_xyt != 1) && (NPE_xzt != 1) && (NPE_yzt != 1)) {
114  vout.crucial(m_vl, "Error at %s: FFTW supports parallelization only in 1 direction of x,y,z.\n",
115  class_name.c_str());
116  exit(EXIT_FAILURE);
117  }
118 
119  bool is_allocated = false;
120  int Lsize_Nsize_prev = 0;
121 
122 #ifdef USE_OPENMP
123  int threads_ok = fftw_init_threads();
124 #endif
125 #ifdef USE_MPI
126  fftw_mpi_init();
127 #endif
128 
129 
130  //- xyz loop
131  for (int i_dim = 0; i_dim < Ndim - 1; ++i_dim) {
132  //- allocate m_in,out = m_in,out[Lsize]
133  if (Lsize_Nsize_prev != Lsize[i_dim] * Nsize[i_dim]) {
134  if (is_allocated) {
135  fftw_free(m_in);
136  fftw_free(m_out);
137  fftw_destroy_plan(m_plan);
138  }
139 
140  is_allocated = false;
141  }
142  Lsize_Nsize_prev = Lsize[i_dim] * Nsize[i_dim];
143 
144  if (!is_allocated) {
145  if (Lsize[i_dim] == Nsize[i_dim]) {
146  size_t fftw_size = sizeof(fftw_complex) * Lsize[i_dim];
147  m_in = (fftw_complex *)fftw_malloc(fftw_size);
148  m_out = (fftw_complex *)fftw_malloc(fftw_size);
149 
150  if (!m_in || !m_out) {
151  vout.crucial(m_vl, "Error at %s: failed to allocate memory %d [Byte].\n",
152  class_name.c_str(), (int)fftw_size);
153  exit(EXIT_FAILURE);
154  }
155  } else {
156 #ifdef USE_MPI
157  ptrdiff_t Lsize_p = Lsize[i_dim];
158  ptrdiff_t fftw_size_p;
159 
160  if (m_is_forward) {
161  fftw_size_p = fftw_mpi_local_size_1d(Lsize_p,
162  Communicator_impl::world(),
163  FFTW_FORWARD, FFTW_ESTIMATE,
164  &m_Nsize_in_p, &m_start_in_p,
165  &m_Nsize_out_p, &m_start_out_p);
166  } else {
167  fftw_size_p = fftw_mpi_local_size_1d(Lsize_p,
168  Communicator_impl::world(),
169  FFTW_BACKWARD, FFTW_ESTIMATE,
170  &m_Nsize_in_p, &m_start_in_p,
171  &m_Nsize_out_p, &m_start_out_p);
172  }
173 
174  m_in = fftw_alloc_complex(fftw_size_p);
175  m_out = fftw_alloc_complex(fftw_size_p);
176 
177  if (!m_in || !m_out) {
178  vout.crucial(m_vl, "Error at %s: failed to allocate memory %d [Byte].\n",
179  class_name.c_str(), (int)fftw_size_p);
180  exit(EXIT_FAILURE);
181  }
182 #endif
183  }
184 
185  is_allocated = true;
186  }
187 
188 
189  //- setup FFTW plan
190 #ifdef USE_OPENMP
191  int Nthread = ThreadManager::get_num_threads();
192  fftw_plan_with_nthreads(Nthread);
193 #endif
194 
195  if (Lsize[i_dim] == Nsize[i_dim]) {
196  if (m_is_forward) {
197  m_plan = fftw_plan_dft_1d(Lsize[i_dim],
198  m_in, m_out,
199  FFTW_FORWARD, FFTW_ESTIMATE);
200  } else {
201  m_plan = fftw_plan_dft_1d(Lsize[i_dim],
202  m_in, m_out,
203  FFTW_BACKWARD, FFTW_ESTIMATE);
204  }
205  } else {
206 #ifdef USE_MPI
207  ptrdiff_t Lsize_p = Lsize[i_dim];
208 
209  if (m_is_forward) {
210  m_plan = fftw_mpi_plan_dft_1d(Lsize_p,
211  m_in, m_out,
212  Communicator_impl::world(),
213  FFTW_FORWARD, FFTW_ESTIMATE);
214  } else {
215  m_plan = fftw_mpi_plan_dft_1d(Lsize_p,
216  m_in, m_out,
217  Communicator_impl::world(),
218  FFTW_BACKWARD, FFTW_ESTIMATE);
219  }
220 #endif
221  }
222 
223 
224  // #### Execution main part ####
225  //- Nin is devided by 2, because of complex(i.e. real and imag)
226  for (int in2 = 0; in2 < Nin / 2; ++in2) {
227  for (int t_global = 0; t_global < Lt; t_global++) {
228  for (int ex = 0; ex < Nex; ++ex) {
229  if (i_dim == 0) {
230  for (int z = 0; z < Nsize[2]; z++) {
231  for (int y = 0; y < Nsize[1]; y++) {
232  for (int xyz_local = 0; xyz_local < Nsize[i_dim]; xyz_local++) {
233  int isite = m_index.site(xyz_local, y, z, t_global);
234  int i_real = 2 * in2;
235  int i_imag = 2 * in2 + 1;
236 
237  m_in[xyz_local][0] = field.cmp(i_real, isite, ex);
238  m_in[xyz_local][1] = field.cmp(i_imag, isite, ex);
239  }
240 
241 
242  fftw_execute(m_plan);
243 
244 
245  for (int xyz_local = 0; xyz_local < Nsize[i_dim]; xyz_local++) {
246  int isite = m_index.site(xyz_local, y, z, t_global);
247  int i_real = 2 * in2;
248  int i_imag = 2 * in2 + 1;
249 
250  field.set(i_real, isite, ex, m_out[xyz_local][0]);
251  field.set(i_imag, isite, ex, m_out[xyz_local][1]);
252  }
253  }
254  }
255  } else if (i_dim == 1) {
256  for (int x = 0; x < Nsize[0]; x++) {
257  for (int z = 0; z < Nsize[2]; z++) {
258  for (int xyz_local = 0; xyz_local < Nsize[i_dim]; xyz_local++) {
259  int isite = m_index.site(x, xyz_local, z, t_global);
260  int i_real = 2 * in2;
261  int i_imag = 2 * in2 + 1;
262 
263  m_in[xyz_local][0] = field.cmp(i_real, isite, ex);
264  m_in[xyz_local][1] = field.cmp(i_imag, isite, ex);
265  }
266 
267 
268  fftw_execute(m_plan);
269 
270 
271  for (int xyz_local = 0; xyz_local < Nsize[i_dim]; xyz_local++) {
272  int isite = m_index.site(x, xyz_local, z, t_global);
273  int i_real = 2 * in2;
274  int i_imag = 2 * in2 + 1;
275 
276  field.set(i_real, isite, ex, m_out[xyz_local][0]);
277  field.set(i_imag, isite, ex, m_out[xyz_local][1]);
278  }
279  }
280  }
281  } else if (i_dim == 2) {
282  for (int y = 0; y < Nsize[1]; y++) {
283  for (int x = 0; x < Nsize[0]; x++) {
284  for (int xyz_local = 0; xyz_local < Nsize[i_dim]; xyz_local++) {
285  int isite = m_index.site(x, y, xyz_local, t_global);
286  int i_real = 2 * in2;
287  int i_imag = 2 * in2 + 1;
288 
289  m_in[xyz_local][0] = field.cmp(i_real, isite, ex);
290  m_in[xyz_local][1] = field.cmp(i_imag, isite, ex);
291  }
292 
293 
294  fftw_execute(m_plan);
295 
296 
297  for (int xyz_local = 0; xyz_local < Nsize[i_dim]; xyz_local++) {
298  int isite = m_index.site(x, y, xyz_local, t_global);
299  int i_real = 2 * in2;
300  int i_imag = 2 * in2 + 1;
301 
302  field.set(i_real, isite, ex, m_out[xyz_local][0]);
303  field.set(i_imag, isite, ex, m_out[xyz_local][1]);
304  }
305  }
306  }
307  }
308  //- end of if( i_dim == 0 ){
309  }
310  }
311  }
312  //- end of outer loops
313  }
314  //- end of for (int i_dim = 0; i_dim < Ndim - 1; ++i_dim)
315 
316  //- normailzation for FFTW_BACKWARD
317  if (!m_is_forward) {
318  scal(field, 1.0 / Lxyz);
319  }
320 
321  //- tidy up
322  if (is_allocated) {
323  fftw_free(m_in);
324  fftw_free(m_out);
325  fftw_destroy_plan(m_plan);
326  }
327 }
328 
329 
330 //====================================================================
331 void FFT_xyz_1dim::fft(Field& field_out, const Field& field_in)
332 {
333  const int Ndim = CommonParameters::Ndim();
334 
335  //- global lattice size
336  std::vector<int> Lsize(Ndim - 1);
337 
338  Lsize[0] = CommonParameters::Lx();
339  Lsize[1] = CommonParameters::Ly();
340  Lsize[2] = CommonParameters::Lz();
341 
342  const int Lt = CommonParameters::Lt();
343  const int Lxyz = Lsize[0] * Lsize[1] * Lsize[2];
344 
345  //- local lattice size
346  std::vector<int> Nsize(Ndim - 1);
347  Nsize[0] = CommonParameters::Nx();
348  Nsize[1] = CommonParameters::Ny();
349  Nsize[2] = CommonParameters::Nz();
350 
351  const int Nin = field_in.nin();
352  const int Nex = field_in.nex();
353 
354  std::vector<int> NPE(Ndim);
355  NPE[0] = CommonParameters::NPEx();
356  NPE[1] = CommonParameters::NPEy();
357  NPE[2] = CommonParameters::NPEz();
358  NPE[3] = CommonParameters::NPEt();
359 
360  const int NPE_xyt = NPE[0] * NPE[1] * NPE[3];
361  const int NPE_xzt = NPE[0] * NPE[2] * NPE[3];
362  const int NPE_yzt = NPE[1] * NPE[2] * NPE[3];
363 
364  if ((NPE_xyt != 1) && (NPE_xzt != 1) && (NPE_yzt != 1)) {
365  vout.crucial(m_vl, "Error at %s: FFTW supports parallelization only in 1 direction.\n",
366  class_name.c_str());
367  exit(EXIT_FAILURE);
368  }
369 
370  bool is_allocated = false;
371  int Lsize_Nsize_prev = 0;
372 
373  field_out = field_in;
374 
375 #ifdef USE_OPENMP
376  int threads_ok = fftw_init_threads();
377 #endif
378 #ifdef USE_MPI
379  fftw_mpi_init();
380 #endif
381 
382 
383  //- xyz loop
384  for (int i_dim = 0; i_dim < Ndim - 1; ++i_dim) {
385  //- allocate m_in,out = m_in,out[Lsize]
386  if (Lsize_Nsize_prev != Lsize[i_dim] * Nsize[i_dim]) {
387  if (is_allocated) {
388  fftw_free(m_in);
389  fftw_free(m_out);
390  fftw_destroy_plan(m_plan);
391  }
392 
393  is_allocated = false;
394  }
395  Lsize_Nsize_prev = Lsize[i_dim] * Nsize[i_dim];
396 
397  if (!is_allocated) {
398  if (Lsize[i_dim] == Nsize[i_dim]) {
399  size_t fftw_size = sizeof(fftw_complex) * Lsize[i_dim];
400  m_in = (fftw_complex *)fftw_malloc(fftw_size);
401  m_out = (fftw_complex *)fftw_malloc(fftw_size);
402 
403  if (!m_in || !m_out) {
404  vout.crucial(m_vl, "Error at %s: failed to allocate memory %d [Byte].\n",
405  class_name.c_str(), (int)fftw_size);
406  exit(EXIT_FAILURE);
407  }
408  } else {
409 #ifdef USE_MPI
410  ptrdiff_t Lsize_p = Lsize[i_dim];
411  ptrdiff_t fftw_size_p;
412 
413  if (m_is_forward) {
414  fftw_size_p = fftw_mpi_local_size_1d(Lsize_p,
415  Communicator_impl::world(),
416  FFTW_FORWARD, FFTW_ESTIMATE,
417  &m_Nsize_in_p, &m_start_in_p,
418  &m_Nsize_out_p, &m_start_out_p);
419  } else {
420  fftw_size_p = fftw_mpi_local_size_1d(Lsize_p,
421  Communicator_impl::world(),
422  FFTW_BACKWARD, FFTW_ESTIMATE,
423  &m_Nsize_in_p, &m_start_in_p,
424  &m_Nsize_out_p, &m_start_out_p);
425  }
426 
427  m_in = fftw_alloc_complex(fftw_size_p);
428  m_out = fftw_alloc_complex(fftw_size_p);
429 
430  if (!m_in || !m_out) {
431  vout.crucial(m_vl, "Error at %s: failed to allocate memory %d [Byte].\n",
432  class_name.c_str(), (int)fftw_size_p);
433  exit(EXIT_FAILURE);
434  }
435 #endif
436  }
437 
438  is_allocated = true;
439  }
440 
441 
442  //- setup FFTW plan
443 #ifdef USE_OPENMP
444  int Nthread = ThreadManager::get_num_threads();
445  fftw_plan_with_nthreads(Nthread);
446 #endif
447 
448  if (Lsize[i_dim] == Nsize[i_dim]) {
449  if (m_is_forward) {
450  m_plan = fftw_plan_dft_1d(Lsize[i_dim],
451  m_in, m_out,
452  FFTW_FORWARD, FFTW_ESTIMATE);
453  } else {
454  m_plan = fftw_plan_dft_1d(Lsize[i_dim],
455  m_in, m_out,
456  FFTW_BACKWARD, FFTW_ESTIMATE);
457  }
458  } else {
459 #ifdef USE_MPI
460  ptrdiff_t Lsize_p = Lsize[i_dim];
461 
462  if (m_is_forward) {
463  m_plan = fftw_mpi_plan_dft_1d(Lsize_p,
464  m_in, m_out,
465  Communicator_impl::world(),
466  FFTW_FORWARD, FFTW_ESTIMATE);
467  } else {
468  m_plan = fftw_mpi_plan_dft_1d(Lsize_p,
469  m_in, m_out,
470  Communicator_impl::world(),
471  FFTW_BACKWARD, FFTW_ESTIMATE);
472  }
473 #endif
474  }
475 
476 
477  // #### Execution main part ####
478  //- Nin is devided by 2, because of complex(i.e. real and imag)
479  for (int in2 = 0; in2 < Nin / 2; ++in2) {
480  for (int t_global = 0; t_global < Lt; t_global++) {
481  for (int ex = 0; ex < Nex; ++ex) {
482  if (i_dim == 0) {
483  for (int z = 0; z < Nsize[2]; z++) {
484  for (int y = 0; y < Nsize[1]; y++) {
485  for (int xyz_local = 0; xyz_local < Nsize[i_dim]; xyz_local++) {
486  int isite = m_index.site(xyz_local, y, z, t_global);
487  int i_real = 2 * in2;
488  int i_imag = 2 * in2 + 1;
489 
490  m_in[xyz_local][0] = field_out.cmp(i_real, isite, ex);
491  m_in[xyz_local][1] = field_out.cmp(i_imag, isite, ex);
492  }
493 
494 
495  fftw_execute(m_plan);
496 
497 
498  for (int xyz_local = 0; xyz_local < Nsize[i_dim]; xyz_local++) {
499  int isite = m_index.site(xyz_local, y, z, t_global);
500  int i_real = 2 * in2;
501  int i_imag = 2 * in2 + 1;
502 
503  field_out.set(i_real, isite, ex, m_out[xyz_local][0]);
504  field_out.set(i_imag, isite, ex, m_out[xyz_local][1]);
505  }
506  }
507  }
508  } else if (i_dim == 1) {
509  for (int x = 0; x < Nsize[0]; x++) {
510  for (int z = 0; z < Nsize[2]; z++) {
511  for (int xyz_local = 0; xyz_local < Nsize[i_dim]; xyz_local++) {
512  int isite = m_index.site(x, xyz_local, z, t_global);
513  int i_real = 2 * in2;
514  int i_imag = 2 * in2 + 1;
515 
516  m_in[xyz_local][0] = field_out.cmp(i_real, isite, ex);
517  m_in[xyz_local][1] = field_out.cmp(i_imag, isite, ex);
518  }
519 
520 
521  fftw_execute(m_plan);
522 
523 
524  for (int xyz_local = 0; xyz_local < Nsize[i_dim]; xyz_local++) {
525  int isite = m_index.site(x, xyz_local, z, t_global);
526  int i_real = 2 * in2;
527  int i_imag = 2 * in2 + 1;
528 
529  field_out.set(i_real, isite, ex, m_out[xyz_local][0]);
530  field_out.set(i_imag, isite, ex, m_out[xyz_local][1]);
531  }
532  }
533  }
534  } else if (i_dim == 2) {
535  for (int y = 0; y < Nsize[1]; y++) {
536  for (int x = 0; x < Nsize[0]; x++) {
537  for (int xyz_local = 0; xyz_local < Nsize[i_dim]; xyz_local++) {
538  int isite = m_index.site(x, y, xyz_local, t_global);
539  int i_real = 2 * in2;
540  int i_imag = 2 * in2 + 1;
541 
542  m_in[xyz_local][0] = field_out.cmp(i_real, isite, ex);
543  m_in[xyz_local][1] = field_out.cmp(i_imag, isite, ex);
544  }
545 
546 
547  fftw_execute(m_plan);
548 
549 
550  for (int xyz_local = 0; xyz_local < Nsize[i_dim]; xyz_local++) {
551  int isite = m_index.site(x, y, xyz_local, t_global);
552  int i_real = 2 * in2;
553  int i_imag = 2 * in2 + 1;
554 
555  field_out.set(i_real, isite, ex, m_out[xyz_local][0]);
556  field_out.set(i_imag, isite, ex, m_out[xyz_local][1]);
557  }
558  }
559  }
560  }
561  //- end of if( i_dim == 0 ){
562  }
563  }
564  }
565  //- end of outer loops
566  }
567  //- end of for (int i_dim = 0; i_dim < Ndim - 1; ++i_dim)
568 
569  //- normailzation for FFTW_BACKWARD
570  if (!m_is_forward) {
571  scal(field_out, 1.0 / Lxyz);
572  }
573 
574  //- tidy up
575  if (is_allocated) {
576  fftw_free(m_in);
577  fftw_free(m_out);
578  fftw_destroy_plan(m_plan);
579  }
580 }
581 
582 
583 //====================================================================
584 void FFT_xyz_1dim::fft(Field& field_out, const Field& field_in, const Direction dir)
585 {
586  // save state
587  bool backup_fwbw = m_is_forward;
588 
589  // find direction and set
590  if (dir == FORWARD) {
591  m_is_forward = true;
592  } else if (dir == BACKWARD) {
593  m_is_forward = false;
594  } else {
595  vout.crucial(m_vl, "%s: unknown direction %d. failed.\n", class_name.c_str(), dir);
596  exit(EXIT_FAILURE);
597  }
598 
599  // delegate to another method
600  fft(field_out, field_in);
601 
602  // restore state
603  m_is_forward = backup_fwbw;
604 }
605 
606 
607 //==========================================================
608 //==================================================END=====
609 #endif
CommonParameters::Ny
static int Ny()
Definition: commonParameters.h:106
CommonParameters::Nz
static int Nz()
Definition: commonParameters.h:107
Parameters::set_string
void set_string(const string &key, const string &value)
Definition: parameters.cpp:39
ThreadManager::get_num_threads
static int get_num_threads()
returns available number of threads.
Definition: threadManager.cpp:246
CommonParameters::Ndim
static int Ndim()
Definition: commonParameters.h:117
Field::set
void set(const int jin, const int site, const int jex, double v)
Definition: field.h:175
Parameters
Class for parameters.
Definition: parameters.h:46
Field::nex
int nex() const
Definition: field.h:128
CommonParameters::Ly
static int Ly()
Definition: commonParameters.h:92
Direction
Direction
Definition: bridge_defs.h:24
Field::nin
int nin() const
Definition: field.h:126
CommonParameters::Nx
static int Nx()
Definition: commonParameters.h:105
CommonParameters::Lt
static int Lt()
Definition: commonParameters.h:94
CommonParameters::Lx
static int Lx()
Definition: commonParameters.h:91
CommonParameters::Lz
static int Lz()
Definition: commonParameters.h:93
CommonParameters::NPEz
static int NPEz()
Definition: commonParameters.h:99
fft_xyz_1dim.h
CommonParameters::NPEy
static int NPEy()
Definition: commonParameters.h:98
Field::cmp
double cmp(const int jin, const int site, const int jex) const
Definition: field.h:143
Bridge::BridgeIO::set_verbose_level
static VerboseLevel set_verbose_level(const std::string &str)
Definition: bridgeIO.cpp:133
CommonParameters::NPEx
static int NPEx()
Definition: commonParameters.h:97
scal
void scal(Field &x, const double a)
scal(x, a): x = a * x
Definition: field.cpp:261
Parameters::fetch_string
int fetch_string(const string &key, string &value) const
Definition: parameters.cpp:378
CommonParameters::NPEt
static int NPEt()
Definition: commonParameters.h:100
Bridge::BridgeIO::crucial
void crucial(const char *format,...)
Definition: bridgeIO.cpp:180
Field
Container of Field-type object.
Definition: field.h:46
Bridge::BridgeIO::general
void general(const char *format,...)
Definition: bridgeIO.cpp:200
Bridge::vout
BridgeIO vout
Definition: bridgeIO.cpp:512
Bridge::BridgeIO::get_verbose_level
static std::string get_verbose_level(const VerboseLevel vl)
Definition: bridgeIO.cpp:154