Functions
QGpCoreWave/FortranInterface.cpp File Reference
#include <QGpCoreTools.h>
#include "FortranInterface.h"
#include "Seismic1DModel.h"
#include "Rayleigh.h"
#include "Dispersion.h"

Functions

void dispersion_curve_init_ (int *verbose)
void dispersion_curve_love_ (int *nLayers, double *h, double *vp, double *vs, double *rho, int *nSamples, double *omega, int *nModes, double *slowness, int *group)
void dispersion_curve_rayleigh_ (int *nLayers, double *h, double *vp, double *vs, double *rho, int *nSamples, double *omega, int *nModes, double *slowness, int *group)

Function Documentation

void dispersion_curve_init_ ( int *  verbose)

If verbose is null, warnings are switched off.

References QGpCoreTools::CoreApplicationPrivate::freezeStream().

{
  new PluginCoreApplication;
  if(*verbose==0) {
    PluginCoreApplication::instance()->freezeStream(true);
  }
}
void dispersion_curve_love_ ( int *  nLayers,
double *  h,
double *  vp,
double *  vs,
double *  rho,
int *  nSamples,
double *  omega,
int *  nModes,
double *  slowness,
int *  group 
)

Computes Love dispersion curve

TODO: Vp is useless... check for a way to avoid it in arguments

References QGpCoreWave::Dispersion::calculate(), QGpCoreWave::Seismic1DModel::initCalculation(), QGpCoreTools::Value< numberType >::isValid(), mode, QGpCoreWave::ModalStorage::mode(), QGpCoreWave::Dispersion::setGroupSlowness(), QGpCoreWave::Seismic1DModel::setH(), QGpCoreWave::Seismic1DModel::setQp(), QGpCoreWave::Seismic1DModel::setQs(), QGpCoreWave::Seismic1DModel::setRho(), QGpCoreWave::Seismic1DModel::setSlowP(), QGpCoreWave::Seismic1DModel::setSlowS(), and QGpCoreTools::Value< numberType >::value().

{
  QGpCoreWave::Seismic1DModel m(*nLayers);
  int nLayers1=*nLayers - 1;
  for(int i=0;i<nLayers1;i++) {
    m.setH(i, h[i]);
    m.setSlowP(i, 1.0/vp[i] );
    m.setSlowS(i, 1.0/vs[i] );
    m.setRho(i, rho[i]);
    m.setQp(i, 0.0);
    m.setQs(i, 0.0);
  }
  m.setSlowP(nLayers1, 1.0/vp[nLayers1] );
  m.setSlowS(nLayers1, 1.0/vs[nLayers1] );
  m.setRho(nLayers1, rho[nLayers1]);
  m.setQp(nLayers1, 0.0);
  m.setQs(nLayers1, 0.0);
  //QTextStream sOut(stdout);
  //m.toStream(sOut);
  m.initCalculation();

  QGpCoreWave::Love love(&m);
  QVector<double> x(*nSamples);
  for(int i=0;i<*nSamples;i++) {
    x[i]=omega[i];
  }
  QGpCoreWave::Dispersion dispersion (*nModes, &x);
  dispersion.calculate(&love, 0);
  if(*group!=0) {
    dispersion.setGroupSlowness();
  }
  for(int iMode=0;iMode<*nModes;iMode++) {
    const RealValue * mode=dispersion.mode(iMode);
    for(int iSample=0; iSample<*nSamples; iSample++) {
      const RealValue& v=mode[iSample];
      if(v.isValid()) {
        *(slowness++)=v.value();
      } else {
        *(slowness++)=0.0;
      }
    }
  }
}
void dispersion_curve_rayleigh_ ( int *  nLayers,
double *  h,
double *  vp,
double *  vs,
double *  rho,
int *  nSamples,
double *  omega,
int *  nModes,
double *  slowness,
int *  group 
)

Computes Rayleigh dispersion curve

References QGpCoreWave::Dispersion::calculate(), QGpCoreWave::Seismic1DModel::initCalculation(), QGpCoreTools::Value< numberType >::isValid(), mode, QGpCoreWave::ModalStorage::mode(), QGpCoreWave::Dispersion::setGroupSlowness(), QGpCoreWave::Seismic1DModel::setH(), QGpCoreWave::Seismic1DModel::setQp(), QGpCoreWave::Seismic1DModel::setQs(), QGpCoreWave::Seismic1DModel::setRho(), QGpCoreWave::Seismic1DModel::setSlowP(), QGpCoreWave::Seismic1DModel::setSlowS(), and QGpCoreTools::Value< numberType >::value().

{
  QGpCoreWave::Seismic1DModel m(*nLayers);
  int nLayers1=*nLayers - 1;
  for(int i=0;i<nLayers1;i++) {
    m.setH(i, h[i]);
    m.setSlowP(i, 1.0/vp[i] );
    m.setSlowS(i, 1.0/vs[i] );
    m.setRho(i, rho[i]);
    m.setQp(i, 0.0);
    m.setQs(i, 0.0);
  }
  m.setSlowP(nLayers1, 1.0/vp[nLayers1] );
  m.setSlowS(nLayers1, 1.0/vs[nLayers1] );
  m.setRho(nLayers1, rho[nLayers1]);
  m.setQp(nLayers1, 0.0);
  m.setQs(nLayers1, 0.0);
  //QTextStream sOut(stdout);
  //m.toStream(sOut);
  m.initCalculation();

  QGpCoreWave::Rayleigh rayleigh(&m);
  //App::stream() << "nSamples=" << *nSamples << ", nModes=" << *nModes << ", group=" << *group << endl;
  QVector<double> x(*nSamples);
  for(int i=0;i<*nSamples;i++) {
    x[i]=omega[i];
    //App::stream() << "frequency[" << i << "]=" << x[i]/(2*M_PI) << " Hz" << endl;
  }
  QGpCoreWave::Dispersion dispersion (*nModes, &x);
  dispersion.calculate(&rayleigh, 0);
  if(*group!=0) {
    dispersion.setGroupSlowness();
  }
  for(int iMode=0;iMode<*nModes;iMode++) {
    //App::stream() << "Mode=" << iMode << endl;
    const RealValue * mode=dispersion.mode(iMode);
    for(int iSample=0; iSample<*nSamples; iSample++) {
      const RealValue& v=mode[iSample];
      if(v.isValid()) {
        *(slowness++)=v.value();
        /*App::stream() << "frequency[" << iSample << "]=" << x[iSample]/(2*M_PI) << " Hz"
                      << ", slowness[" << iSample << "]=" << v.value()
                      << ", velocity[" << iSample << "]=" << 1.0/v.value() << endl;*/
      } else {
        *(slowness++)=0.0;
      }
    }
  }
}
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Defines