Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / applications / utilities / postProcessing / patch / patchIntegrate / patchIntegrate.C
blob4c55b18aa5b1124db0a55a67005c823c163444b3
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  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"
33 #include "cyclicPolyPatch.H"
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 // Main program:
38 int main(int argc, char *argv[])
40 #   include "addRegionOption.H"
41     timeSelector::addOptions();
42     argList::validArgs.append("fieldName");
43     argList::validArgs.append("patchName");
44 #   include "setRootCase.H"
45 #   include "createTime.H"
46     instantList timeDirs = timeSelector::select0(runTime, args);
47 #   include "createNamedMesh.H"
49     word fieldName(args.additionalArgs()[0]);
50     word patchName(args.additionalArgs()[1]);
52     forAll(timeDirs, timeI)
53     {
54         runTime.setTime(timeDirs[timeI], timeI);
55         Info<< "Time = " << runTime.timeName() << endl;
57         IOobject fieldHeader
58         (
59             fieldName,
60             runTime.timeName(),
61             mesh,
62             IOobject::MUST_READ
63         );
65         // Check field exists
66         if (fieldHeader.headerOk())
67         {
68             mesh.readUpdate();
70             label patchi = mesh.boundaryMesh().findPatchID(patchName);
71             if (patchi < 0)
72             {
73                 FatalError
74                     << "Unable to find patch " << patchName << nl
75                     << exit(FatalError);
76             }
78             // Give patch area
79             if (isA<cyclicPolyPatch>(mesh.boundaryMesh()[patchi]))
80             {
81                 Info<< "    Cyclic patch vector area: " << nl;
82                 label nFaces = mesh.boundaryMesh()[patchi].size();
83                 vector sum1 = vector::zero;
84                 vector sum2 = vector::zero;
85                 for (label i=0; i<nFaces/2; i++)
86                 {
87                     sum1 += mesh.Sf().boundaryField()[patchi][i];
88                     sum2 += mesh.Sf().boundaryField()[patchi][i+nFaces/2];
89                 }
90                 reduce(sum1, sumOp<vector>());
91                 reduce(sum2, sumOp<vector>());
92                 Info<< "    - half 1 = " << sum1 << ", " << mag(sum1) << nl
93                     << "    - half 2 = " << sum2 << ", " << mag(sum2) << nl
94                     << "    - total  = " << (sum1 + sum2) << ", "
95                     << mag(sum1 + sum2) << endl;
96                 Info<< "    Cyclic patch area magnitude = "
97                     << gSum(mesh.magSf().boundaryField()[patchi])/2.0 << endl;
98             }
99             else
100             {
101                 Info<< "    Area vector of patch "
102                     << patchName << '[' << patchi << ']' << " = "
103                     << gSum(mesh.Sf().boundaryField()[patchi]) << endl;
104                 Info<< "    Area magnitude of patch "
105                     << patchName << '[' << patchi << ']' << " = "
106                     << gSum(mesh.magSf().boundaryField()[patchi]) << endl;
107             }
109             // Read field and calc integral
110             if (fieldHeader.headerClassName() == volScalarField::typeName)
111             {
112                 Info<< "    Reading " << volScalarField::typeName << " "
113                     << fieldName << endl;
115                 volScalarField field(fieldHeader, mesh);
117                 Info<< "    Integral of " << fieldName
118                     << " over vector area of patch "
119                     << patchName << '[' << patchi << ']' << " = "
120                     << gSum
121                        (
122                            mesh.Sf().boundaryField()[patchi]
123                           *field.boundaryField()[patchi]
124                        )
125                     << nl;
127                 Info<< "    Integral of " << fieldName
128                     << " over area magnitude of patch "
129                     << patchName << '[' << patchi << ']' << " = "
130                     << gSum
131                        (
132                            mesh.magSf().boundaryField()[patchi]
133                           *field.boundaryField()[patchi]
134                        )
135                     << nl;
136             }
137             else if
138             (
139                 fieldHeader.headerClassName() == surfaceScalarField::typeName
140             )
141             {
142                 Info<< "    Reading " << surfaceScalarField::typeName << " "
143                     << fieldName << endl;
145                 surfaceScalarField field(fieldHeader, mesh);
146                 scalar sumField = gSum(field.boundaryField()[patchi]);
148                 Info<< "    Integral of " << fieldName << " over patch "
149                     << patchName << '[' << patchi << ']' << " = "
150                     << sumField << nl;
151             }
152             else
153             {
154                 FatalError
155                     << "Only possible to integrate "
156                     << volScalarField::typeName << "s "
157                     << "and " << surfaceScalarField::typeName << "s"
158                     << nl << exit(FatalError);
159             }
160         }
161         else
162         {
163             Info<< "    No field " << fieldName << endl;
164         }
166         Info<< endl;
167     }
169     Info<< "End\n" << endl;
171     return 0;
175 // ************************************************************************* //