Changeset 1405 in Sophya for trunk


Ignore:
Timestamp:
Feb 15, 2001, 6:58:16 PM (25 years ago)
Author:
ansari
Message:

Ajout documentation - Reza 15/2/2001

Location:
trunk
Files:
10 edited

Legend:

Unmodified
Added
Removed
  • trunk/SophyaExt/IFFTW/fftwserver.cc

    r1403 r1405  
    33#include "FFTW/rfftw.h"
    44
     5/*!
     6  \class SOPHYA::FFTWServer
     7  \ingroup IFFTW
     8  An implementation of FFTServerInterface based on FFTW, double
     9  precision arrays, using FFTW package, availabale from http://www.fftw.org.
     10  Refer to FFTServerInterface for details about FFTServer operations.
     11
     12  \code
     13  #include "fftwserver.h"
     14  // ...
     15  TMatrix<r_8> in(24,32);
     16  TMatrix< complex<r_8> > out;
     17  in = RandomSequence();
     18  FFTWServer ffts;
     19  ffts.setNormalize(true);  // To have normalized transforms
     20  cout << " FFTServer info string= " << ffts.getInfo() << endl;
     21  cout << "in= " << in << endl;
     22  cout << " Calling ffts.FFTForward(in, out) : " << endl;
     23  ffts.FFTForward(in, out);
     24  cout << "out= " << out << endl;
     25  \endcode
     26
     27*/
    528
    629#define MAXND_FFTW 5
     30
     31namespace SOPHYA {
    732
    833class FFTWServerPlan {
     
    2752   
    2853};
     54
     55} // Fin du namespace
    2956
    3057FFTWServerPlan::FFTWServerPlan(int n, fftw_direction dir, bool fgreal)
     
    4067  _n = n; 
    4168  _dir = dir;
    42   if (fgreal) rp = rfftw_create_plan(n, dir, FFTW_ESTIMATE);
    43   else p = fftw_create_plan(n, dir, FFTW_ESTIMATE);
     69  if (fgreal) {
     70    rp = rfftw_create_plan(n, dir, FFTW_ESTIMATE);
     71    if (rp == NULL)
     72      throw AllocationError("FFTWServerPlan: failed to create plan (rp) !");
     73  }
     74  else {
     75    p = fftw_create_plan(n, dir, FFTW_ESTIMATE);
     76    if (p == NULL)
     77      throw AllocationError("FFTWServerPlan: failed to create plan (p) !");
     78  }
     79   
    4480}
    4581
     
    6298  for(k=nd; k<MAXND_FFTW; k++) _nxyz[k] = -10;
    6399  _dir = dir;
    64   if (fgreal) rpnd = rfftwnd_create_plan(_nd, _nxyz, dir, FFTW_ESTIMATE);
    65   else pnd = fftwnd_create_plan(_nd, _nxyz, dir, FFTW_ESTIMATE);
     100  if (fgreal) {
     101    rpnd = rfftwnd_create_plan(_nd, _nxyz, dir, FFTW_ESTIMATE);
     102    if (rpnd == NULL)
     103      throw AllocationError("FFTWServerPlan: failed to create plan (rpnd) !");
     104  }
     105  else {
     106    pnd = fftwnd_create_plan(_nd, _nxyz, dir, FFTW_ESTIMATE);
     107    if (pnd == NULL)
     108      throw AllocationError("FFTWServerPlan: failed to create plan (pnd) !");
     109  }
     110
    66111}
    67112
     
    86131    fftw_destroy_plan(p);
    87132    p = fftw_create_plan(n, _dir, FFTW_ESTIMATE);
     133    if (p == NULL)
     134      throw AllocationError("FFTWServerPlan::Recreate failed to create plan (p) !");
    88135  }
    89136  else {
    90137    rfftw_destroy_plan(rp);
    91138    rp = rfftw_create_plan(n, _dir, FFTW_ESTIMATE);
    92   }
     139    if (rp == NULL)
     140      throw AllocationError("FFTWServerPlan::Recreate failed to create plan (rp) !");
     141  }
     142
    93143}
    94144
     
    117167    fftwnd_destroy_plan(pnd);
    118168    pnd = fftwnd_create_plan(_nd, _nxyz, _dir, FFTW_ESTIMATE);
     169    if (pnd == NULL)
     170      throw AllocationError("FFTWServerPlan::Recreate failed to create plan (pnd) !");
    119171  }
    120172  else {
    121173    rfftwnd_destroy_plan(rpnd);
    122174    rpnd = rfftwnd_create_plan(_nd, _nxyz, _dir, FFTW_ESTIMATE);
     175    if (rpnd == NULL)
     176      throw AllocationError("FFTWServerPlan::Recreate failed to create plan (rpnd) !");
    123177  }
    124178
  • trunk/SophyaExt/IFFTW/fftwserver.h

    r1403 r1405  
    66
    77
    8 // Classe definissant l'interface pour les transformees de Fourier
    9 // L'implementation par defaut est vide et lance une exception
     8// Classe implementant l'interface FFTServerInterface en
     9// utilisant FFTW
     10
     11namespace SOPHYA {
    1012
    1113class FFTWServerPlan;
     
    4446};
    4547
     48} // Fin du namespace
     49
    4650#endif
  • trunk/SophyaLib/BaseTools/dvlist.cc

    r1386 r1405  
    5353   used to represent FITS headers. The class \b SOPHYA::ObjFileIO<DVList>
    5454   handles serialisation for DVList. (See SOPHYA::PPersist ).
    55    
     55   \sa SOPHYA::ObjFileIO<DVList>
     56
     57   \code
    5658   //  ------- Using DVList objects ------
    5759   DVList  dvl;
     
    7678
    7779   \endcode
    78    \sa SOPHYA::ObjFileIO<DVList>
    7980*/
    8081
  • trunk/SophyaLib/HiStats/ntuple.cc

    r1371 r1405  
    1616   Simple NTuple class (a Table or 2-D data set) with floating value
    1717   columns.
     18   \sa SOPHYA::ObjFileIO<NTuple>
     19
     20   \code
     21   #include "ntuple.h"
     22   // ...
     23   char * names[3] = {"XPos", "YPos", "Val"};
     24   // NTuple (Table) creation with 3 columns
     25   NTuple  nt(3, names);
     26   // Filling the NTuple
     27   r_4 x[3];
     28   for(int i=0; i<63; i++) {
     29     x[0] = (i%9)-4.;  x[1] = (i/9)-3.;  x[2] = x[0]*x[0]+x[1]*x[1];
     30     nt.Fill(x);
     31   }
     32   // Printing table info
     33   cout << nt ;
     34   // Saving object into a PPF file
     35   POutPersist po("nt.ppf");
     36   po << nt ;
     37   \endcode
    1838*/
    1939
     
    3252
    3353/* --Methode-- */
     54//! Default constructor
    3455//++
    3556NTuple::NTuple()
     
    4667
    4768
     69//! Constructor with specification of number of columns and column name
     70/*!
     71  \param nvar : Number of columns in the table
     72  \param noms : Array of column names (length(name) < 8 characters)
     73  \param blk : Optional argument specifying number of table lines
     74  in a given data block
     75 */
    4876//++
    4977NTuple::NTuple(int nvar, char** noms, int blk)
     
    77105}
    78106
     107//! Copy constructor - Copies the table definition and associated data
    79108//                                       cmv 8/10/99
    80109//++
     
    108137
    109138/* --Methode-- */
     139//! Constructor with table initialized from a PPF file
    110140//++
    111141NTuple::NTuple(char *flnm)
     
    131161
    132162/* --Methode-- */
     163//! Clear the data table definition and deletes the associated data
    133164void NTuple::Clean()
    134165{
     
    149180
    150181/* --Methode--        cmv 08/10/99 */
     182//! = operator, copies the data table definition and its contents
    151183//++
    152184NTuple& NTuple::operator = (const NTuple& NT)
     
    181213
    182214/* --Methode-- */
     215//! Adds an entry (a line) to the table
     216/*!
     217  \param x : content of the line to be appended to the table
     218 */
    183219//++
    184220void  NTuple::Fill(r_4* x)
     
    317353
    318354/* --Methode-- */
     355//! Prints table definition and number of entries
    319356//++
    320357void  NTuple::Show(ostream& os)  const
     
    339376
    340377/* --Methode-- */
     378//! Fills the table, reading lines from an ascii file
     379/*!
     380  \param fn : file name
     381  \param defval : default value for empty cells in the ascii file
     382 */
    341383//++
    342384int  NTuple::FillFromASCIIFile(string const& fn, float defval)
     
    544586}
    545587
     588/*!
     589   \class SOPHYA::ObjFileIO<NTuple>
     590   \ingroup HiStats
     591   Persistence (serialisation) handler for class NTuple
     592*/
    546593
    547594/* --Methode-- */
  • trunk/SophyaLib/HiStats/ntuple.h

    r1046 r1405  
    9797};
    9898
     99/*! Prints table information on stream \b s (nt.Show(os)) */
    99100inline ostream& operator << (ostream& s, NTuple const & nt)
    100101  {  nt.Show(s);  return(s);  }
    101102
     103/*! Writes the object in the POutPersist stream \b os */
    102104inline POutPersist& operator << (POutPersist& os, NTuple & obj)
    103105{ ObjFileIO<NTuple> fio(&obj);  fio.Write(os);  return(os); }
     106/*! Reads the object from the PInPersist stream \b is */
    104107inline PInPersist& operator >> (PInPersist& is, NTuple & obj)
    105108{ ObjFileIO<NTuple> fio(&obj);  fio.Read(is);  return(is); }
  • trunk/SophyaLib/HiStats/xntuple.cc

    r1371 r1405  
    2020   columns with int, float or string type content. In addition, this class
    2121   can handle large data sets, using swap space on disk.
     22   \sa SOPHYA::ObjFileIO<XNTuple>
     23
     24   \code
     25   #include "xntuple.h"
     26   // ...
     27   char * names[4] = {"X", "X2", "XInt","XStr"};
     28   // XNTuple (Table) creation with 4 columns, of integer,
     29   // double(2) and string type
     30   XNTuple  xnt(2,0,1,1, names);
     31   // Filling the NTuple
     32   r_8 xd[2];
     33   int_4 xi[2];
     34   char xss[2][32];
     35   char * xs[2] = {xss[0], xss[1]} ;
     36   for(int i=0; i<50; i++) {
     37     xi[0] = i;  xd[0] = i+0.5;  xd[1] = xd[0]*xd[0];
     38     sprintf(xs[0],"X=%g", xd[0]);
     39     xnt.Fill(xd, NULL, xi, xs);
     40   }
     41   // Printing table info
     42   cout << xnt ;
     43   // Saving object into a PPF file
     44   POutPersist po("xnt.ppf");
     45   po << xnt ;
     46   \endcode
     47   
    2248*/
    2349
     
    218244//--
    219245
     246//! Constructor with specification colum types and names
     247/*!
     248  \param ndvar : Number of columns with type double (r_8)
     249  \param nfvar : Number of columns with type float (r_4)
     250  \param nivar : Number of columns with type integer (int_4)
     251  \param nsvar : Number of columns with type string
     252  \param vnames : Column names (in the order r_8 r_4 int_4 string)
     253  (colum names are limited to 31 characters)
     254  \param blk : data block size (in number of table lines)
     255  \param maxblk : Maximum number of data block in memory
     256  \param strsz : string length (for string type columns)
     257 */
    220258XNTuple::XNTuple(int ndvar, int nfvar, int nivar, int nsvar,
    221259                   char** vnames,
     
    249287
    250288
     289//! Constructor with table initialized from a PPF file
    251290XNTuple::XNTuple(string const& flnm)
    252291    : mNEnt(0), mNBlk(0), mNSwBlk(0),
     
    268307
    269308
     309//! Copy constructor - Copies the table definition and associated data
    270310XNTuple::XNTuple(XNTuple const& nt)
    271311    : mNEnt(0), mNBlk(0), mNSwBlk(0),
     
    294334
    295335
     336//! Clear the data table definition and deletes the associated data
    296337void XNTuple::clean()
    297338{
     
    341382//      Opérateur égal (=) , copie "nt" dans le premier NTuple
    342383//--
     384//! Appends an entry (line) to the table
     385/*!
     386  \param d_data : double (r_8) line elements
     387  \param f_data : float (r_4) line elements
     388  \param i_data : integer (int_4) line elements
     389  \param s_data : string line elements
     390 */
    343391void XNTuple::Fill(r_8* d_data, r_4* f_data, int_4* i_data, char** s_data)
    344392{
     
    552600//--
    553601
     602//! Returns the associated DVList object
    554603DVList&  XNTuple::Info()
    555604{
     
    597646
    598647
     648//! Prints table definition and number of entries
    599649void   XNTuple::Show(ostream& os) const
    600650{
     
    636686}
    637687
     688//! NOT YET IMPLEMENTED ! - Fills table from an ascii file
    638689int    XNTuple::FillFromASCIIFile(string const& fn, r_8 ddval,  r_4 dfval,
    639690                                        int dival, const char * dsval)
     
    646697
    647698
     699//! Defines swap file path for XNTuple objects (Default=/tmp)
    648700void  XNTuple::SetSwapPath(char* p)
    649701{
     
    843895}
    844896
     897/*!
     898   \class SOPHYA::ObjFileIO<XNTuple>
     899   \ingroup HiStats
     900   Persistence (serialisation) handler for class XNTuple
     901*/
    845902
    846903void   ObjFileIO<XNTuple>::WriteSelf(POutPersist& ppout) const
  • trunk/SophyaLib/HiStats/xntuple.h

    r1046 r1405  
    159159};
    160160
    161 
     161/*! Prints table information on stream \b s (nt.Show(os)) */
    162162inline ostream& operator << (ostream& s, XNTuple const & nt)
    163163  {  nt.Show(s);  return(s);  }
    164164
     165/*! Writes the object in the POutPersist stream \b os */
    165166inline POutPersist& operator << (POutPersist& os, XNTuple & obj)
    166167{ ObjFileIO<XNTuple> fio(&obj);  fio.Write(os);  return(os); }
     168/*! Reads the object from the PInPersist stream \b is */
    167169inline PInPersist& operator >> (PInPersist& is, XNTuple & obj)
    168170{ ObjFileIO<XNTuple> fio(&obj);  fio.Read(is);  return(is); }
  • trunk/SophyaLib/NTools/fftpserver.cc

    r1402 r1405  
    1515necessarily correspond with the equivalent fftpack function.  For example,
    1616fftpack "forward" transformations are in fact inverse fourier transformations.
    17 Otherwise, the output is in the fftpack format.
    18 
    1917
    2018Due to the way that fftpack manages
     
    2321of differing length, it may be more efficient to create an fftserver object
    2422for each length.
     23
     24  \code
     25  #include "fftpserver.h"
     26  // ...
     27  TVector<r_8> in(32);
     28  TVector< complex<r_8> > out;
     29  in = RandomSequence();
     30  FFTPackServer ffts;
     31  ffts.setNormalize(true);  // To have normalized transforms
     32  cout << " FFTServer info string= " << ffts.getInfo() << endl;
     33  cout << "in= " << in << endl;
     34  cout << " Calling ffts.FFTForward(in, out) : " << endl;
     35  ffts.FFTForward(in, out);
     36  cout << "out= " << out << endl;
     37  \endcode
    2538*/
    2639
  • trunk/SophyaLib/NTools/fftservintf.cc

    r1402 r1405  
    66  \ingroup NTools
    77  Defines the interface for FFT (Fast Fourier Transform) operations.
     8  Definitions :
     9    - Sampling period \b T
     10    - Sampling frequency \b fs=1/T
     11    - Total number of samples \b N
     12    - Frequency step in Fourier space \b =fs/N=1/(N*T)
     13    - Component frequencies
     14        - k=0      ->  0
     15        - k=1      ->  1/(N*T)
     16        - k        ->  k/(N*T)
     17        - k=N/2    ->  1/(2*T)   (Nyquist frequency)
     18        - k>N/2    ->  k/(N*T)   (or negative frequency -(N-k)/(N*T))
     19
     20  For a sampling period T=1, the computed Fourier components correspond to :
     21  \verbatim
     22  0  1/N  2/N  ... 1/2  1/2+1/N  1/2+2/N ... 1-2/N  1-1/N
     23  0  1/N  2/N  ... 1/2                   ...  -2/N   -1/N
     24  \endverbatim
     25
     26  For complex one-dimensional transforms:
     27  \f[
     28  out(i) = F_{norm} \Sigma_{j} \ e^{-2 \pi \sqrt{-1} \ i \  j} \ {\rm (forward)}
     29  \f]
     30  \f[
     31  out(i) = F_{norm} \Sigma_{j} \ e^{2 \pi \sqrt{-1} \ i \  j} \ {\rm (backward)}
     32  \f]
     33  i,j= 0..N-1 , where N is the input or the output array size.
     34
     35  For complex multi-dimensional transforms:
     36  \f[
     37  out(i1,i2,...,id) = F_{norm} \Sigma_{j1} \Sigma_{j2} ... \Sigma_{jd} \
     38  e^{-2 \pi \sqrt{-1} \ i1 \ j1} ... e^{-2 \pi \sqrt{-1} \ id \ jd} \ {\rm (forward)}
     39  \f]
     40  \f[
     41  out(i1,i2,...,id) = F_{norm} \Sigma_{j1} \Sigma_{j2} ... \Sigma_{jd} \
     42  e^{2 \pi \sqrt{-1} \ i1 \ j1} ... e^{2 \pi \sqrt{-1} \ id \ jd} \ {\rm (backward)}
     43  \f]
     44
     45  For real forward transforms, the input array is real, and
     46  the output array complex, with Fourier components up to k=N/2.
     47  For real backward transforms, the input array is complex and
     48  the output array is real.
    849*/
    950
     
    2364
    2465/* --Methode-- */
     66//! Forward Fourier transform for double precision complex data
     67/*!
     68  \param in : Input complex array
     69  \param out : Output complex array
     70 */
    2571void FFTServerInterface::FFTForward(TArray< complex<r_8> > const &, TArray< complex<r_8> > &)
    2672{
     
    2975
    3076/* --Methode-- */
     77//! Backward (inverse) Fourier transform for double precision complex data
     78/*!
     79  \param in : Input complex array
     80  \param out : Output complex array
     81 */
    3182void FFTServerInterface::FFTBackward(TArray< complex<r_8> > const &, TArray< complex<r_8> > &)
    3283{
     
    3586
    3687/* --Methode-- */
     88//! Forward Fourier transform for double precision real input data
     89/*!
     90  \param in : Input real array
     91  \param out : Output complex array
     92 */
    3793void FFTServerInterface::FFTForward(TArray< r_8 > const &, TArray< complex<r_8> > &)
    3894{
     
    4197
    4298/* --Methode-- */
     99//! Backward (inverse) Fourier transform for double precision real output data
     100/*!
     101  \param in : Input complex array
     102  \param out : Output real array
     103  \param usoutsz : if true, use the output array size for computing the inverse FFT.
     104 */
    43105void FFTServerInterface::FFTBackward(TArray< complex<r_8> > const &, TArray< r_8 > &, bool)
    44106{
     
    50112
    51113/* --Methode-- */
     114//! Forward Fourier transform for complex data
     115/*!
     116  \param in : Input complex array
     117  \param out : Output complex array
     118 */
    52119void FFTServerInterface::FFTForward(TArray< complex<r_4> > const &, TArray< complex<r_4> > &)
    53120{
     
    56123
    57124/* --Methode-- */
     125//! Backward (inverse) Fourier transform for complex data
     126/*!
     127  \param in : Input complex array
     128  \param out : Output complex array
     129 */
    58130void FFTServerInterface::FFTBackward(TArray< complex<r_4> > const &, TArray< complex<r_4> > &)
    59131{
     
    62134
    63135/* --Methode-- */
     136//! Forward Fourier transform for real input data
     137/*!
     138  \param in : Input real array
     139  \param out : Output complex array
     140 */
    64141void FFTServerInterface::FFTForward(TArray< r_4 > const &, TArray< complex<r_4> > &)
    65142{
     
    68145
    69146/* --Methode-- */
     147//! Backward (inverse) Fourier transform for real output data
     148/*!
     149  \param in : Input complex array
     150  \param out : Output real array
     151  \param usoutsz : if true, use the output array size for computing the inverse FFT.
     152 */
    70153void FFTServerInterface::FFTBackward(TArray< complex<r_4> > const &, TArray< r_4 > &, bool)
    71154{
     
    76159
    77160/* --Methode-- */
     161/*!
     162  \class SOPHYA::FFTArrayChecker
     163  \ingroup NTools
     164  Service class for checking array size and resizing output arrays,
     165  to be used by FFTServer classes
     166*/
     167
    78168template <class T>
    79169FFTArrayChecker<T>::FFTArrayChecker(string msg, bool checkpack, bool onedonly)
  • trunk/SophyaLib/NTools/fftservintf.h

    r1402 r1405  
    2727  virtual FFTServerInterface * Clone() = 0;
    2828
     29//! Set/clear the flag for normalizing Fourier transforms.
    2930  inline void setNormalize(bool fg=false) { _fgnorm = fg; }
     31//! Returns the status of normalization flag for the server
    3032  inline bool getNormalize() const { return(_fgnorm); }
     33//! Returns the information string associated with the server
    3134  inline string getInfo() const { return _info; }
    3235
     
    5356};
    5457
    55 } // Fin du namespace
    5658
    5759template <class T>
     
    7981};
    8082
     83} // Fin du namespace
     84
    8185#endif
Note: See TracChangeset for help on using the changeset viewer.