File indexing completed on 2025-01-05 04:12:50
0001 /*************************************************************************** 0002 * * 0003 * copyright : (C) 2007 The University of Toronto * 0004 * netterfield@astro.utoronto.ca * 0005 * copyright : (C) 2005 University of British Columbia * 0006 * dscott@phas.ubc.ca * 0007 * * 0008 * This program is free software; you can redistribute it and/or modify * 0009 * it under the terms of the GNU General Public License as published by * 0010 * the Free Software Foundation; either version 2 of the License, or * 0011 * (at your option) any later version. * 0012 * * 0013 ***************************************************************************/ 0014 0015 0016 #include "autocorrelation.h" 0017 #include "objectstore.h" 0018 #include "ui_autocorrelationconfig.h" 0019 0020 #include <gsl/gsl_fft_real.h> 0021 #include <gsl/gsl_fft_halfcomplex.h> 0022 0023 static const QString& VECTOR_IN = "Vector In"; 0024 static const QString& VECTOR_OUT_AUTO = "Auto-Correlated"; 0025 static const QString& VECTOR_OUT_STEP = "Step Value"; 0026 0027 class ConfigAutoCorrelationPlugin : public Kst::DataObjectConfigWidget, public Ui_AutoCorrelationConfig { 0028 public: 0029 ConfigAutoCorrelationPlugin(QSettings* cfg) : DataObjectConfigWidget(cfg), Ui_AutoCorrelationConfig() { 0030 _store = 0; 0031 setupUi(this); 0032 } 0033 0034 ~ConfigAutoCorrelationPlugin() {} 0035 0036 void setObjectStore(Kst::ObjectStore* store) { 0037 _store = store; 0038 _vector->setObjectStore(store); 0039 } 0040 0041 void setupSlots(QWidget* dialog) { 0042 if (dialog) { 0043 connect(_vector, SIGNAL(selectionChanged(QString)), dialog, SIGNAL(modified())); 0044 } 0045 } 0046 0047 Kst::VectorPtr selectedVector() { return _vector->selectedVector(); }; 0048 void setSelectedVector(Kst::VectorPtr vector) { return _vector->setSelectedVector(vector); }; 0049 0050 virtual void setupFromObject(Kst::Object* dataObject) { 0051 if (AutoCorrelationSource* source = static_cast<AutoCorrelationSource*>(dataObject)) { 0052 setSelectedVector(source->vector()); 0053 } 0054 } 0055 0056 virtual bool configurePropertiesFromXml(Kst::ObjectStore *store, QXmlStreamAttributes& attrs) { 0057 Q_UNUSED(store); 0058 Q_UNUSED(attrs); 0059 0060 bool validTag = true; 0061 0062 // QStringRef av; 0063 // av = attrs.value("value"); 0064 // if (!av.isNull()) { 0065 // _configValue = QVariant(av.toString()).toBool(); 0066 // } 0067 0068 return validTag; 0069 } 0070 0071 public slots: 0072 virtual void save() { 0073 if (_cfg) { 0074 _cfg->beginGroup("Auto Correlation DataObject Plugin"); 0075 _cfg->setValue("Input Vector", _vector->selectedVector()->Name()); 0076 _cfg->endGroup(); 0077 } 0078 } 0079 0080 virtual void load() { 0081 if (_cfg && _store) { 0082 _cfg->beginGroup("Auto Correlation DataObject Plugin"); 0083 QString vectorName = _cfg->value("Input Vector").toString(); 0084 Kst::Object* object = _store->retrieveObject(vectorName); 0085 Kst::Vector* vector = static_cast<Kst::Vector*>(object); 0086 if (vector) { 0087 setSelectedVector(vector); 0088 } 0089 _cfg->endGroup(); 0090 } 0091 } 0092 0093 private: 0094 Kst::ObjectStore *_store; 0095 0096 }; 0097 0098 0099 AutoCorrelationSource::AutoCorrelationSource(Kst::ObjectStore *store) 0100 : Kst::BasicPlugin(store) { 0101 } 0102 0103 0104 AutoCorrelationSource::~AutoCorrelationSource() { 0105 } 0106 0107 0108 QString AutoCorrelationSource::_automaticDescriptiveName() const { 0109 return tr("Auto Correlation Plugin Object"); 0110 } 0111 0112 0113 void AutoCorrelationSource::change(Kst::DataObjectConfigWidget *configWidget) { 0114 if (ConfigAutoCorrelationPlugin* config = static_cast<ConfigAutoCorrelationPlugin*>(configWidget)) { 0115 setInputVector(VECTOR_IN, config->selectedVector()); 0116 } 0117 } 0118 0119 0120 void AutoCorrelationSource::setupOutputs() { 0121 setOutputVector(VECTOR_OUT_AUTO, ""); 0122 setOutputVector(VECTOR_OUT_STEP, ""); 0123 } 0124 0125 0126 bool AutoCorrelationSource::algorithm() { 0127 Kst::VectorPtr inputVector = _inputVectors[VECTOR_IN]; 0128 Kst::VectorPtr outputVectorAuto = _outputVectors[VECTOR_OUT_AUTO]; 0129 Kst::VectorPtr outputVectorStep = _outputVectors[VECTOR_OUT_STEP]; 0130 0131 //Make sure there is at least 1 element in the input vector 0132 if (inputVector->length() < 1) { 0133 _errorString = tr("Error: Input Vector invalid size"); 0134 return false; 0135 } 0136 0137 double* pdArrayOne; 0138 double* pdResult; 0139 double* pdCorrelate; 0140 double dReal; 0141 double dImag; 0142 0143 int iLength; 0144 int iLengthNew; 0145 0146 bool bReturn = false; 0147 0148 // 0149 // zero-pad the array... 0150 // 0151 iLength = inputVector->length(); 0152 iLength *= 2; 0153 0154 outputVectorAuto->resize(inputVector->length(), false); 0155 outputVectorStep->resize(inputVector->length(), false); 0156 0157 // 0158 // round iLength up to the nearest power of two... 0159 // 0160 iLengthNew = 64; 0161 while( iLengthNew < iLength && iLengthNew > 0) { 0162 iLengthNew *= 2; 0163 } 0164 iLength = iLengthNew; 0165 0166 if (iLength <= 0) { 0167 _errorString = tr("Error: Invalid Input length calculated"); 0168 return false; 0169 } 0170 0171 pdArrayOne = new double[iLength]; 0172 if (pdArrayOne != NULL) { 0173 // 0174 // zero-pad the two arrays... 0175 // 0176 memset( pdArrayOne, 0, iLength * sizeof( double ) ); 0177 memcpy( pdArrayOne, inputVector->noNanValue(), inputVector->length() * sizeof( double ) ); 0178 0179 // 0180 // calculate the FFTs of the two functions... 0181 // 0182 if (gsl_fft_real_radix2_transform( pdArrayOne, 1, iLength ) == 0) { 0183 // 0184 // multiply the FFT by its complex conjugate... 0185 // 0186 for (int i=0; i<iLength/2; i++) { 0187 if (i==0 || i==(iLength/2)-1) { 0188 pdArrayOne[i] *= pdArrayOne[i]; 0189 } else { 0190 dReal = pdArrayOne[i] * pdArrayOne[i] + pdArrayOne[iLength-i] * pdArrayOne[iLength-i]; 0191 dImag = pdArrayOne[i] * pdArrayOne[iLength-i] - pdArrayOne[iLength-i] * pdArrayOne[i]; 0192 0193 pdArrayOne[i] = dReal; 0194 pdArrayOne[iLength-i] = dImag; 0195 } 0196 } 0197 0198 // 0199 // do the inverse FFT... 0200 // 0201 if (gsl_fft_halfcomplex_radix2_inverse( pdArrayOne, 1, iLength ) == 0) { 0202 pdResult = outputVectorStep->raw_V_ptr(); 0203 pdCorrelate = outputVectorAuto->raw_V_ptr(); 0204 0205 if (pdResult != NULL && pdCorrelate != NULL) { 0206 for (int i = 0; i < inputVector->length(); ++i) { 0207 outputVectorStep->raw_V_ptr()[i] = pdResult[i]; 0208 } 0209 for (int i = 0; i < inputVector->length(); ++i) { 0210 outputVectorAuto->raw_V_ptr()[i] = pdCorrelate[i]; 0211 } 0212 0213 for (int i = 0; i < inputVector->length(); i++) { 0214 outputVectorStep->raw_V_ptr()[i] = (double)( i - ( inputVector->length() / 2 ) ); 0215 } 0216 0217 memcpy( &(outputVectorAuto->raw_V_ptr()[inputVector->length() / 2]), 0218 &(pdArrayOne[0]), 0219 ( ( inputVector->length() + 1 ) / 2 ) * sizeof( double ) ); 0220 0221 memcpy( &(outputVectorAuto->raw_V_ptr()[0]), 0222 &(pdArrayOne[iLength - (inputVector->length() / 2)]), 0223 ( inputVector->length() / 2 ) * sizeof( double ) ); 0224 0225 bReturn = true; 0226 } 0227 } 0228 } 0229 } 0230 delete[] pdArrayOne; 0231 0232 return bReturn; 0233 } 0234 0235 0236 Kst::VectorPtr AutoCorrelationSource::vector() const { 0237 return _inputVectors[VECTOR_IN]; 0238 } 0239 0240 0241 QStringList AutoCorrelationSource::inputVectorList() const { 0242 return QStringList( VECTOR_IN ); 0243 } 0244 0245 0246 QStringList AutoCorrelationSource::inputScalarList() const { 0247 return QStringList( /*SCALAR_IN*/ ); 0248 } 0249 0250 0251 QStringList AutoCorrelationSource::inputStringList() const { 0252 return QStringList( /*STRING_IN*/ ); 0253 } 0254 0255 0256 QStringList AutoCorrelationSource::outputVectorList() const { 0257 QStringList vectors(VECTOR_OUT_AUTO); 0258 vectors += VECTOR_OUT_STEP; 0259 return vectors; 0260 } 0261 0262 0263 QStringList AutoCorrelationSource::outputScalarList() const { 0264 return QStringList( /*SCALAR_OUT*/ ); 0265 } 0266 0267 0268 QStringList AutoCorrelationSource::outputStringList() const { 0269 return QStringList( /*STRING_OUT*/ ); 0270 } 0271 0272 0273 void AutoCorrelationSource::saveProperties(QXmlStreamWriter &s) { 0274 Q_UNUSED(s); 0275 // s.writeAttribute("value", _configValue); 0276 } 0277 0278 0279 QString AutoCorrelationPlugin::pluginName() const { return tr("Auto Correlation"); } 0280 QString AutoCorrelationPlugin::pluginDescription() const { return tr("Generates the auto-correlation of a vector."); } 0281 0282 0283 Kst::DataObject *AutoCorrelationPlugin::create(Kst::ObjectStore *store, Kst::DataObjectConfigWidget *configWidget, bool setupInputsOutputs) const { 0284 0285 if (ConfigAutoCorrelationPlugin* config = static_cast<ConfigAutoCorrelationPlugin*>(configWidget)) { 0286 0287 AutoCorrelationSource* object = store->createObject<AutoCorrelationSource>(); 0288 0289 if (setupInputsOutputs) { 0290 object->setupOutputs(); 0291 object->setInputVector(VECTOR_IN, config->selectedVector()); 0292 } 0293 0294 object->setPluginName(pluginName()); 0295 0296 object->writeLock(); 0297 object->registerChange(); 0298 object->unlock(); 0299 0300 return object; 0301 } 0302 return 0; 0303 } 0304 0305 0306 Kst::DataObjectConfigWidget *AutoCorrelationPlugin::configWidget(QSettings *settingsObject) const { 0307 ConfigAutoCorrelationPlugin *widget = new ConfigAutoCorrelationPlugin(settingsObject); 0308 return widget; 0309 } 0310 0311 #ifndef QT5 0312 Q_EXPORT_PLUGIN2(kstplugin_BinPlugin, AutoCorrelationPlugin) 0313 #endif 0314 0315 // vim: ts=2 sw=2 et