Page MenuHomec4science

statistics.cpp
No OneTemporary

File Metadata

Created
Thu, May 9, 23:04

statistics.cpp

/**
* @file
* @section LICENSE
*
* Copyright (©) 2016-19 EPFL (École Polytechnique Fédérale de Lausanne),
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Affero General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU Affero General Public License for more details.
*
* You should have received a copy of the GNU Affero General Public License
* along with this program. If not, see <https://www.gnu.org/licenses/>.
*
*/
/* -------------------------------------------------------------------------- */
#include "statistics.hh"
#include "loop.hh"
#include "static_types.hh"
/* -------------------------------------------------------------------------- */
namespace tamaas {
template <UInt dim>
Real Statistics<dim>::computeRMSHeights(Grid<Real, dim>& surface) {
return std::sqrt(surface.var());
}
template <UInt dim>
Real Statistics<dim>::computeSpectralRMSSlope(Grid<Real, dim>& surface) {
auto h_size = GridHermitian<Real, dim>::hermitianDimensions(surface.sizes());
auto wavevectors =
FFTransform<Real, dim>::template computeFrequencies<true>(h_size);
wavevectors *= 2 * M_PI; // need q for slopes
auto psd = computePowerSpectrum(surface);
Real rms_slope_mean = Loop::reduce<operation::plus>(
[] CUDA_LAMBDA(auto q, const auto& psd_val) {
// Checking if we're in the zone that does not have hermitian symmetry
if (std::abs(q.back()) < 1e-15)
return q.l2squared() * psd_val.real();
else
return 2 * q.l2squared() * psd_val.real();
},
range<VectorProxy<Real, dim>>(wavevectors), psd);
return std::sqrt(rms_slope_mean);
}
/* -------------------------------------------------------------------------- */
template <UInt dim>
GridHermitian<Real, dim>
Statistics<dim>::computePowerSpectrum(Grid<Real, dim>& surface) {
auto h_size = GridHermitian<Real, dim>::hermitianDimensions(surface.sizes());
GridHermitian<Real, dim> psd(h_size, surface.getNbComponents());
FFTPlanManager::get().createPlan(surface, psd).forwardTransform();
Real factor = 1. / surface.getNbPoints();
// Squaring the fourier transform of surface and normalizing
Loop::loop(
[factor] CUDA_LAMBDA(Complex & c) {
c *= factor;
c *= conj(c);
},
psd);
FFTPlanManager::get().destroyPlan(surface, psd);
return psd;
}
/* -------------------------------------------------------------------------- */
template <UInt dim>
Grid<Real, dim>
Statistics<dim>::computeAutocorrelation(Grid<Real, dim>& surface) {
Grid<Real, dim> acf(surface.sizes(), surface.getNbComponents());
auto psd = computePowerSpectrum(surface);
FFTPlanManager::get().createPlan(acf, psd).backwardTransform();
acf *= acf.getNbPoints();
FFTPlanManager::get().destroyPlan(acf, psd);
return acf;
}
/* -------------------------------------------------------------------------- */
namespace {
template <UInt dim>
class moment_helper {
public:
moment_helper(const std::array<UInt, dim>& exp) : exponent(exp) {}
CUDA_LAMBDA Complex operator()(VectorProxy<Real, dim> q,
const Complex& phi) const {
Real mul = 1;
for (UInt i = 0; i < dim; ++i)
mul *= std::pow(q(i), exponent[i]);
// Do not duplicate everything from hermitian symmetry
if (std::abs(q.back()) < 1e-15)
return mul * phi;
else
return 2 * mul * phi;
}
private:
std::array<UInt, dim> exponent;
};
} // namespace
template <>
std::vector<Real> Statistics<1>::computeMoments(Grid<Real, 1>& surface) {
constexpr UInt dim = 1;
std::vector<Real> moments(3);
const auto psd = computePowerSpectrum(surface);
auto wavevectors =
FFTransform<Real, dim>::template computeFrequencies<true>(psd.sizes());
// we don't multiply by 2 pi because moments are computed with k
moments[0] = Loop::reduce<operation::plus>(moment_helper<dim>{{{0}}},
range<PVector>(wavevectors), psd)
.real();
moments[1] = Loop::reduce<operation::plus>(moment_helper<dim>{{{2}}},
range<PVector>(wavevectors), psd)
.real();
moments[2] = Loop::reduce<operation::plus>(moment_helper<dim>{{{4}}},
range<PVector>(wavevectors), psd)
.real();
return moments;
}
template <>
std::vector<Real> Statistics<2>::computeMoments(Grid<Real, 2>& surface) {
constexpr UInt dim = 2;
std::vector<Real> moments(3);
const auto psd = computePowerSpectrum(surface);
auto wavevectors =
FFTransform<Real, dim>::template computeFrequencies<true>(psd.sizes());
// we don't multiply by 2 pi because moments are computed with k
moments[0] = Loop::reduce<operation::plus>(moment_helper<dim>{{{0, 0}}},
range<PVector>(wavevectors), psd)
.real();
auto m02 = Loop::reduce<operation::plus>(moment_helper<dim>{{{0, 2}}},
range<PVector>(wavevectors), psd)
.real();
auto m20 = Loop::reduce<operation::plus>(moment_helper<dim>{{{2, 0}}},
range<PVector>(wavevectors), psd)
.real();
moments[1] = 0.5 * (m02 + m20);
auto m22 = Loop::reduce<operation::plus>(moment_helper<dim>{{{2, 2}}},
range<PVector>(wavevectors), psd)
.real();
auto m40 = Loop::reduce<operation::plus>(moment_helper<dim>{{{4, 0}}},
range<PVector>(wavevectors), psd)
.real();
auto m04 = Loop::reduce<operation::plus>(moment_helper<dim>{{{0, 4}}},
range<PVector>(wavevectors), psd)
.real();
moments[2] = (3 * m22 + m40 + m04) / 3.;
return moments;
}
template struct Statistics<1>;
template struct Statistics<2>;
} // namespace tamaas

Event Timeline