BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / applications / utilities / postProcessing / patch / patchIntegrate / patchIntegrate.C
blob8696c6ffde90c9f5a48a65fa29ea6e69b3cf851d
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 Application
25     patchIntegrate
27 Description
28     Calculates the integral of the specified field over the specified patch.
30 \*---------------------------------------------------------------------------*/
32 #include "fvCFD.H"
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 // Main program:
37 int main(int argc, char *argv[])
39 #   include "addRegionOption.H"
40     timeSelector::addOptions();
41     argList::validArgs.append("fieldName");
42     argList::validArgs.append("patchName");
43 #   include "setRootCase.H"
44 #   include "createTime.H"
45     instantList timeDirs = timeSelector::select0(runTime, args);
46 #   include "createNamedMesh.H"
48     const word fieldName = args[1];
49     const word patchName = args[2];
51     forAll(timeDirs, timeI)
52     {
53         runTime.setTime(timeDirs[timeI], timeI);
54         Info<< "Time = " << runTime.timeName() << endl;
56         IOobject fieldHeader
57         (
58             fieldName,
59             runTime.timeName(),
60             mesh,
61             IOobject::MUST_READ
62         );
64         // Check field exists
65         if (fieldHeader.headerOk())
66         {
67             mesh.readUpdate();
69             const label patchI = mesh.boundaryMesh().findPatchID(patchName);
70             if (patchI < 0)
71             {
72                 FatalError
73                     << "Unable to find patch " << patchName << nl
74                     << exit(FatalError);
75             }
77             // Give patch area
78             Info<< "    Area vector of patch "
79                 << patchName << '[' << patchI << ']' << " = "
80                 << gSum(mesh.Sf().boundaryField()[patchI]) << endl;
81             Info<< "    Area magnitude of patch "
82                 << patchName << '[' << patchI << ']' << " = "
83                 << gSum(mesh.magSf().boundaryField()[patchI]) << endl;
85             // Read field and calc integral
86             if (fieldHeader.headerClassName() == volScalarField::typeName)
87             {
88                 Info<< "    Reading " << volScalarField::typeName << " "
89                     << fieldName << endl;
91                 volScalarField field(fieldHeader, mesh);
93                 Info<< "    Integral of " << fieldName
94                     << " over vector area of patch "
95                     << patchName << '[' << patchI << ']' << " = "
96                     << gSum
97                        (
98                            mesh.Sf().boundaryField()[patchI]
99                           *field.boundaryField()[patchI]
100                        )
101                     << nl;
103                 Info<< "    Integral of " << fieldName
104                     << " over area magnitude of patch "
105                     << patchName << '[' << patchI << ']' << " = "
106                     << gSum
107                        (
108                            mesh.magSf().boundaryField()[patchI]
109                           *field.boundaryField()[patchI]
110                        )
111                     << nl;
112             }
113             else if
114             (
115                 fieldHeader.headerClassName() == surfaceScalarField::typeName
116             )
117             {
118                 Info<< "    Reading " << surfaceScalarField::typeName << " "
119                     << fieldName << endl;
121                 surfaceScalarField field(fieldHeader, mesh);
122                 scalar sumField = gSum(field.boundaryField()[patchI]);
124                 Info<< "    Integral of " << fieldName << " over patch "
125                     << patchName << '[' << patchI << ']' << " = "
126                     << sumField << nl;
127             }
128             else
129             {
130                 FatalError
131                     << "Only possible to integrate "
132                     << volScalarField::typeName << "s "
133                     << "and " << surfaceScalarField::typeName << "s"
134                     << nl << exit(FatalError);
135             }
136         }
137         else
138         {
139             Info<< "    No field " << fieldName << endl;
140         }
142         Info<< endl;
143     }
145     Info<< "End\n" << endl;
147     return 0;
151 // ************************************************************************* //