Prevent reference XML reader segfault on nullptr error strings
[gromacs.git] / src / gromacs / correlationfunctions / crosscorr.cpp
blob22291ee132ee0e9bd902458b3a0e89b4be1fc876
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2017, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
35 /*! \internal
36 * \file
37 * \brief
38 * Implements routine for computing a cross correlation between two data sets
40 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
41 * \ingroup module_correlationfunctions
43 #include "gmxpre.h"
45 #include "crosscorr.h"
47 #include "gromacs/fft/fft.h"
48 #include "gromacs/utility/exceptions.h"
49 #include "gromacs/utility/smalloc.h"
51 /*! \brief
52 * return the size witch the array should be after zero padding
54 * \param[in] n Factor to multiply by 2
55 * \return zeroPaddingSize
57 static int zeroPaddingSize(int n)
59 return 2*n;
62 /*! \brief
63 * Compute complex conjugate. Output in the first input variable.
65 * \param[in] in1 first complex number
66 * \param[in] in2 second complex number
68 static void complexConjugatMult(t_complex *in1, t_complex *in2)
70 t_complex res;
71 res.re = in1->re * in2->re + in1->im * in2->im;
72 res.im = in1->re * -in2->im + in1->im * in2->re;
73 in1->re = res.re;
74 in1->im = res.im;
77 /*! \brief
78 * Compute one cross correlation corr = f x g using FFT.
80 * \param[in] n number of data point
81 * \param[in] f first function
82 * \param[in] g second function
83 * \param[out] corr output correlation
84 * \param[in] fft FFT data structure
86 static void cross_corr_low(int n, real f[], real g[], real corr[], gmx_fft_t fft)
88 int i;
89 const int size = zeroPaddingSize(n);
90 t_complex * in1, * in2;
92 snew(in1, size);
93 snew(in2, size);
95 for (i = 0; i < n; i++)
97 in1[i].re = f[i];
98 in1[i].im = 0;
99 in2[i].re = g[i];
100 in2[i].im = 0;
102 for (; i < size; i++)
104 in1[i].re = 0;
105 in1[i].im = 0;
106 in2[i].re = 0;
107 in2[i].im = 0;
109 gmx_fft_1d(fft, GMX_FFT_FORWARD, in1, in1);
110 gmx_fft_1d(fft, GMX_FFT_FORWARD, in2, in2);
112 for (i = 0; i < size; i++)
114 complexConjugatMult(&in1[i], &in2[i]);
115 in1[i].re /= size;
117 gmx_fft_1d(fft, GMX_FFT_BACKWARD, in1, in1);
119 for (i = 0; i < n; i++)
121 corr[i] = in1[i].re;
124 sfree(in1);
125 sfree(in2);
129 void cross_corr(int n, real f[], real g[], real corr[])
131 gmx_fft_t fft;
132 gmx_fft_init_1d(&fft, zeroPaddingSize(n), GMX_FFT_FLAG_CONSERVATIVE);
133 cross_corr_low( n, f, g, corr, fft);
134 gmx_fft_destroy(fft);
135 gmx_fft_cleanup();
138 void many_cross_corr(int nFunc, int * nData, real ** f, real ** g, real ** corr)
140 #pragma omp parallel
141 //gmx_fft_t is not thread safe, so structure are allocated per thread.
143 int i;
145 #pragma omp for
146 for (i = 0; i < nFunc; i++)
150 gmx_fft_t fft;
151 gmx_fft_init_1d(&fft, zeroPaddingSize(nData[i]), GMX_FFT_FLAG_CONSERVATIVE);
152 cross_corr_low( nData[i], f[i], g[i], corr[i], fft);
153 gmx_fft_destroy(fft);
155 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
158 gmx_fft_cleanup();