File: noise.C

package info (click to toggle)
openfoam 4.1%2Bdfsg1-1
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 163,028 kB
  • ctags: 58,990
  • sloc: cpp: 830,760; sh: 10,227; ansic: 8,215; xml: 745; lex: 437; awk: 194; sed: 91; makefile: 77; python: 18
file content (185 lines) | stat: -rw-r--r-- 5,778 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2011-2016 OpenFOAM Foundation
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.

    OpenFOAM 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 General Public License
    for more details.

    You should have received a copy of the GNU General Public License
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

Application
    noise

Description
    Utility to perform noise analysis of pressure data using the noiseFFT
    library.

    Control settings are read from the $FOAM_CASE/system/noiseDict dictionary,
    or user-specified dictionary using the -dict option.  Pressure data is
    read using a CSV reader:

Usage
    \verbatim
    pRef        101325;
    N           65536;
    nw          100;
    f1          25;
    fU          10000;
    graphFormat raw;

    pressureData
    {
        fileName        "pressureData"
        nHeaderLine         1;          // number of header lines
        refColumn           0;          // reference column index
        componentColumns    (1);        // component column indices
        separator           " ";        // optional (defaults to ",")
        mergeSeparators     no;         // merge multiple separators
        outOfBounds         clamp;      // optional out-of-bounds handling
        interpolationScheme linear;     // optional interpolation scheme
    }
    \endverbatim

    where
    \table
        Property    | Description                   | Required  | Default value
        pRef        | Reference pressure            | no        | 0
        N           | Number of samples in sampling window | no | 65536
        nw          | Number of sampling windows    | no        | 100
        fl          | Lower frequency band          | no        | 25
        fU          | Upper frequency band          | no        | 10000
        graphFormat | Output graph format          | no        | raw
    \endtable

    Current graph outputs include:
    - FFT of the pressure data
    - narrow-band PFL (pressure-fluctuation level) spectrum
    - one-third-octave-band PFL spectrum
    - one-third-octave-band pressure spectrum

See also
    CSV.H
    noiseFFT.H

\*---------------------------------------------------------------------------*/


#include "noiseFFT.H"
#include "argList.H"
#include "Time.H"
#include "writeFiles.H"
#include "CSV.H"

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

using namespace Foam;

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

Foam::scalar checkUniformTimeStep(const scalarField& t)
{
    // check that a uniform time step has been applied
    scalar deltaT = -1.0;
    if (t.size() > 1)
    {
        for (label i = 1; i < t.size(); i++)
        {
            scalar dT = t[i] - t[i-1];
            if (deltaT < 0)
            {
                deltaT = dT;
            }

            if (mag(deltaT - dT) > SMALL)
            {
                FatalErrorInFunction
                    << "Unable to process data with a variable time step"
                    << exit(FatalError);
            }
        }
    }
    else
    {
        FatalErrorInFunction
            << "Unable to create FFT with a single value"
            << exit(FatalError);
    }

    return deltaT;
}


// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

int main(int argc, char *argv[])
{
    argList::noParallel();
    #include "addDictOption.H"
    #include "setRootCase.H"
    #include "createTime.H"
    #include "createFields.H"

    Info<< "Reading data file" << endl;
    Function1Types::CSV<scalar> pData("pressure", dict, "Data");

    // time history data
    const scalarField t(pData.x());

    // pressure data
    const scalarField p(pData.y());

    if (t.size() < N)
    {
        FatalErrorInFunction
            << "Block size N = " << N
            << " is larger than number of data = " << t.size()
            << exit(FatalError);
    }

    Info<< "    read " << t.size() << " values" << nl << endl;


    Info<< "Creating noise FFT" << endl;
    noiseFFT nfft(checkUniformTimeStep(t), p);

    nfft -= pRef;

    fileName baseFileName(pData.fName().lessExt());

    graph Pf(nfft.RMSmeanPf(N, min(nfft.size()/N, nw)));
    Info<< "    Creating graph for " << Pf.title() << endl;
    Pf.write(baseFileName + graph::wordify(Pf.title()), graphFormat);

    graph Lf(nfft.Lf(Pf));
    Info<< "    Creating graph for " << Lf.title() << endl;
    Lf.write(baseFileName + graph::wordify(Lf.title()), graphFormat);

    graph Ldelta(nfft.Ldelta(Lf, f1, fU));
    Info<< "    Creating graph for " << Ldelta.title() << endl;
    Ldelta.write(baseFileName + graph::wordify(Ldelta.title()), graphFormat);

    graph Pdelta(nfft.Pdelta(Pf, f1, fU));
    Info<< "    Creating graph for " << Pdelta.title() << endl;
    Pdelta.write(baseFileName + graph::wordify(Pdelta.title()), graphFormat);

    Info<< nl << "End\n" << endl;

    return 0;
}


// ************************************************************************* //