BUGFIX: Uninitialised member variables
[foam-extend-3.2.git] / applications / utilities / postProcessing / patch / patchIntegrate / patchIntegrate.C
blob8d03368b2fd4e5f152701e1a0839622234cc2747
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Application
26     patchIntegrate
28 Description
29     Calculates the integral of the specified field over the specified patch.
31 \*---------------------------------------------------------------------------*/
33 #include "fvCFD.H"
34 #include "cyclicPolyPatch.H"
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 // Main program:
39 int main(int argc, char *argv[])
41 #   include "addRegionOption.H"
42     timeSelector::addOptions();
43     argList::validArgs.append("fieldName");
44     argList::validArgs.append("patchName");
45 #   include "setRootCase.H"
46 #   include "createTime.H"
47     instantList timeDirs = timeSelector::select0(runTime, args);
48 #   include "createNamedMesh.H"
50     word fieldName(args.additionalArgs()[0]);
51     word patchName(args.additionalArgs()[1]);
53     forAll(timeDirs, timeI)
54     {
55         runTime.setTime(timeDirs[timeI], timeI);
56         Info<< "Time = " << runTime.timeName() << endl;
58         IOobject fieldHeader
59         (
60             fieldName,
61             runTime.timeName(),
62             mesh,
63             IOobject::MUST_READ
64         );
66         // Check field exists
67         if (fieldHeader.headerOk())
68         {
69             mesh.readUpdate();
71             label patchi = mesh.boundaryMesh().findPatchID(patchName);
72             if (patchi < 0)
73             {
74                 FatalError
75                     << "Unable to find patch " << patchName << nl
76                     << exit(FatalError);
77             }
79             // Give patch area
80             if (isA<cyclicPolyPatch>(mesh.boundaryMesh()[patchi]))
81             {
82                 Info<< "    Cyclic patch vector area: " << nl;
83                 label nFaces = mesh.boundaryMesh()[patchi].size();
84                 vector sum1 = vector::zero;
85                 vector sum2 = vector::zero;
86                 for (label i=0; i<nFaces/2; i++)
87                 {
88                     sum1 += mesh.Sf().boundaryField()[patchi][i];
89                     sum2 += mesh.Sf().boundaryField()[patchi][i+nFaces/2];
90                 }
91                 reduce(sum1, sumOp<vector>());
92                 reduce(sum2, sumOp<vector>());
93                 Info<< "    - half 1 = " << sum1 << ", " << mag(sum1) << nl
94                     << "    - half 2 = " << sum2 << ", " << mag(sum2) << nl
95                     << "    - total  = " << (sum1 + sum2) << ", "
96                     << mag(sum1 + sum2) << endl;
97                 Info<< "    Cyclic patch area magnitude = "
98                     << gSum(mesh.magSf().boundaryField()[patchi])/2.0 << endl;
99             }
100             else
101             {
102                 Info<< "    Area vector of patch "
103                     << patchName << '[' << patchi << ']' << " = "
104                     << gSum(mesh.Sf().boundaryField()[patchi]) << endl;
105                 Info<< "    Area magnitude of patch "
106                     << patchName << '[' << patchi << ']' << " = "
107                     << gSum(mesh.magSf().boundaryField()[patchi]) << endl;
108             }
110             // Read field and calc integral
111             if (fieldHeader.headerClassName() == volScalarField::typeName)
112             {
113                 Info<< "    Reading " << volScalarField::typeName << " "
114                     << fieldName << endl;
116                 volScalarField field(fieldHeader, mesh);
118                 Info<< "    Integral of " << fieldName
119                     << " over vector area of patch "
120                     << patchName << '[' << patchi << ']' << " = "
121                     << gSum
122                        (
123                            mesh.Sf().boundaryField()[patchi]
124                           *field.boundaryField()[patchi]
125                        )
126                     << nl;
128                 Info<< "    Integral of " << fieldName
129                     << " over area magnitude of patch "
130                     << patchName << '[' << patchi << ']' << " = "
131                     << gSum
132                        (
133                            mesh.magSf().boundaryField()[patchi]
134                           *field.boundaryField()[patchi]
135                        )
136                     << nl;
137             }
138             else if
139             (
140                 fieldHeader.headerClassName() == surfaceScalarField::typeName
141             )
142             {
143                 Info<< "    Reading " << surfaceScalarField::typeName << " "
144                     << fieldName << endl;
146                 surfaceScalarField field(fieldHeader, mesh);
147                 scalar sumField = gSum(field.boundaryField()[patchi]);
149                 Info<< "    Integral of " << fieldName << " over patch "
150                     << patchName << '[' << patchi << ']' << " = "
151                     << sumField << nl;
152             }
153             else
154             {
155                 FatalError
156                     << "Only possible to integrate "
157                     << volScalarField::typeName << "s "
158                     << "and " << surfaceScalarField::typeName << "s"
159                     << nl << exit(FatalError);
160             }
161         }
162         else
163         {
164             Info<< "    No field " << fieldName << endl;
165         }
167         Info<< endl;
168     }
170     Info<< "End\n" << endl;
172     return 0;
176 // ************************************************************************* //