BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / postProcessing / functionObjects / field / fieldValues / faceSource / faceSource.C
blob5141604eed0e8094f75b549f5a02c8495cd1d686
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 \*---------------------------------------------------------------------------*/
26 #include "faceSource.H"
27 #include "fvMesh.H"
28 #include "cyclicPolyPatch.H"
29 #include "emptyPolyPatch.H"
30 #include "coupledPolyPatch.H"
31 #include "sampledSurface.H"
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 defineTypeNameAndDebug(Foam::fieldValues::faceSource, 0);
37 namespace Foam
39     template<>
40     const char* Foam::NamedEnum
41     <
42         Foam::fieldValues::faceSource::sourceType,
43         3
44     >::names[] =
45     {
46         "faceZone",
47         "patch",
48         "sampledSurface"
49     };
52     template<>
53     const char* Foam::NamedEnum
54     <
55         Foam::fieldValues::faceSource::operationType,
56         7
57     >::names[] =
58     {
59         "none",
60         "sum",
61         "areaAverage",
62         "areaIntegrate",
63         "weightedAverage",
64         "min",
65         "max"
66     };
71 const Foam::NamedEnum<Foam::fieldValues::faceSource::sourceType, 3>
72     Foam::fieldValues::faceSource::sourceTypeNames_;
74 const Foam::NamedEnum<Foam::fieldValues::faceSource::operationType, 7>
75     Foam::fieldValues::faceSource::operationTypeNames_;
78 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
80 void Foam::fieldValues::faceSource::setFaceZoneFaces()
82     label zoneId = mesh().faceZones().findZoneID(sourceName_);
84     if (zoneId < 0)
85     {
86         FatalErrorIn("faceSource::faceSource::setFaceZoneFaces()")
87             << type() << " " << name_ << ": "
88             << sourceTypeNames_[source_] << "(" << sourceName_ << "):" << nl
89             << "    Unknown face zone name: " << sourceName_
90             << ". Valid face zones are: " << mesh().faceZones().names()
91             << nl << exit(FatalError);
92     }
94     const faceZone& fZone = mesh().faceZones()[zoneId];
96     faceId_.setSize(fZone.size());
97     facePatchId_.setSize(fZone.size());
98     faceSign_.setSize(fZone.size());
100     label count = 0;
101     forAll(fZone, i)
102     {
103         label faceI = fZone[i];
105         label faceId = -1;
106         label facePatchId = -1;
107         if (mesh().isInternalFace(faceI))
108         {
109             faceId = faceI;
110             facePatchId = -1;
111         }
112         else
113         {
114             facePatchId = mesh().boundaryMesh().whichPatch(faceI);
115             const polyPatch& pp = mesh().boundaryMesh()[facePatchId];
116             if (isA<coupledPolyPatch>(pp))
117             {
118                 if (refCast<const coupledPolyPatch>(pp).owner())
119                 {
120                     faceId = pp.whichFace(faceI);
121                 }
122                 else
123                 {
124                     faceId = -1;
125                 }
126             }
127             else if (!isA<emptyPolyPatch>(pp))
128             {
129                 faceId = faceI - pp.start();
130             }
131             else
132             {
133                 faceId = -1;
134                 facePatchId = -1;
135             }
136         }
138         if (faceId >= 0)
139         {
140             if (fZone.flipMap()[i])
141             {
142                 faceSign_[count] = -1;
143             }
144             else
145             {
146                 faceSign_[count] = 1;
147             }
148             faceId_[count] = faceId;
149             facePatchId_[count] = facePatchId;
150             count++;
151         }
152     }
154     faceId_.setSize(count);
155     facePatchId_.setSize(count);
156     faceSign_.setSize(count);
157     nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
159     if (debug)
160     {
161         Pout<< "Original face zone size = " << fZone.size()
162             << ", new size = " << count << endl;
163     }
167 void Foam::fieldValues::faceSource::setPatchFaces()
169     const label patchId = mesh().boundaryMesh().findPatchID(sourceName_);
171     if (patchId < 0)
172     {
173         FatalErrorIn("faceSource::constructFaceAddressing()")
174             << type() << " " << name_ << ": "
175             << sourceTypeNames_[source_] << "(" << sourceName_ << "):" << nl
176             << "    Unknown patch name: " << sourceName_
177             << ". Valid patch names are: "
178             << mesh().boundaryMesh().names() << nl
179             << exit(FatalError);
180     }
182     const polyPatch& pp = mesh().boundaryMesh()[patchId];
184     label nFaces = pp.size();
185     if (isA<cyclicPolyPatch>(pp))
186     {
187         nFaces /= 2;
188     }
189     else if (isA<emptyPolyPatch>(pp))
190     {
191         nFaces = 0;
192     }
194     faceId_.setSize(nFaces);
195     facePatchId_.setSize(nFaces);
196     faceSign_.setSize(nFaces);
197     nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
199     forAll(faceId_, faceI)
200     {
201         faceId_[faceI] = faceI;
202         facePatchId_[faceI] = patchId;
203         faceSign_[faceI] = 1;
204     }
208 void Foam::fieldValues::faceSource::sampledSurfaceFaces(const dictionary& dict)
210     surfacePtr_ = sampledSurface::New
211     (
212         name_,
213         mesh(),
214         dict.subDict("sampledSurfaceDict")
215     );
216     surfacePtr_().update();
217     nFaces_ = returnReduce(surfacePtr_().faces().size(), sumOp<label>());
221 // * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
223 void Foam::fieldValues::faceSource::initialise(const dictionary& dict)
225     switch (source_)
226     {
227         case stFaceZone:
228         {
229             setFaceZoneFaces();
230             break;
231         }
232         case stPatch:
233         {
234             setPatchFaces();
235             break;
236         }
237         case stSampledSurface:
238         {
239             sampledSurfaceFaces(dict);
240             break;
241         }
242         default:
243         {
244             FatalErrorIn("faceSource::initialise()")
245                 << type() << " " << name_ << ": "
246                 << sourceTypeNames_[source_] << "(" << sourceName_ << "):"
247                 << nl << "    Unknown source type. Valid source types are:"
248                 << sourceTypeNames_.sortedToc() << nl << exit(FatalError);
249         }
250     }
252     if (nFaces_ == 0)
253     {
254         WarningIn
255         (
256             "Foam::fieldValues::faceSource::initialise(const dictionary&)"
257         )
258             << type() << " " << name_ << ": "
259             << sourceTypeNames_[source_] << "(" << sourceName_ << "):" << nl
260             << "    Source has no faces - deactivating" << endl;
262         active_ = false;
263         return;
264     }
266     scalar totalArea;
268     if (surfacePtr_.valid())
269     {
270         surfacePtr_().update();
271         totalArea = gSum(surfacePtr_().magSf());
272     }
273     else
274     {
275         totalArea = gSum(filterField(mesh().magSf(), false));
276     }
278     Info<< type() << " " << name_ << ":" << nl
279         << "    total faces  = " << nFaces_
280         << nl
281         << "    total area   = " << totalArea
282         << nl;
284     if (operation_ == opWeightedAverage)
285     {
286         dict.lookup("weightField") >> weightFieldName_;
287         Info<< "    weight field = " << weightFieldName_;
288     }
290     Info<< nl << endl;
294 void Foam::fieldValues::faceSource::writeFileHeader()
296     if (outputFilePtr_.valid())
297     {
298         outputFilePtr_()
299             << "# Source : " << sourceTypeNames_[source_] << " "
300             << sourceName_ <<  nl << "# Faces  : " << nFaces_ << nl
301             << "# Time" << tab << "sum(magSf)";
303         forAll(fields_, i)
304         {
305             outputFilePtr_()
306                 << tab << operationTypeNames_[operation_]
307                 << "(" << fields_[i] << ")";
308         }
310         outputFilePtr_() << endl;
311     }
315 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
317 Foam::fieldValues::faceSource::faceSource
319     const word& name,
320     const objectRegistry& obr,
321     const dictionary& dict,
322     const bool loadFromFiles
325     fieldValue(name, obr, dict, loadFromFiles),
326     source_(sourceTypeNames_.read(dict.lookup("source"))),
327     operation_(operationTypeNames_.read(dict.lookup("operation"))),
328     weightFieldName_("undefinedWeightedFieldName"),
329     nFaces_(0),
330     faceId_(),
331     facePatchId_(),
332     faceSign_()
334     read(dict);
338 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
340 Foam::fieldValues::faceSource::~faceSource()
344 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
346 void Foam::fieldValues::faceSource::read(const dictionary& dict)
348     fieldValue::read(dict);
350     if (active_)
351     {
352         initialise(dict);
353     }
357 void Foam::fieldValues::faceSource::write()
359     fieldValue::write();
361     if (active_)
362     {
363         scalar totalArea;
365         if (surfacePtr_.valid())
366         {
367             surfacePtr_().update();
368             totalArea = gSum(surfacePtr_().magSf());
369         }
370         else
371         {
372             totalArea = gSum(filterField(mesh().magSf(), false));
373         }
376         if (Pstream::master())
377         {
378             outputFilePtr_() << obr_.time().value() << tab << totalArea;
379         }
381         forAll(fields_, i)
382         {
383             writeValues<scalar>(fields_[i]);
384             writeValues<vector>(fields_[i]);
385             writeValues<sphericalTensor>(fields_[i]);
386             writeValues<symmTensor>(fields_[i]);
387             writeValues<tensor>(fields_[i]);
388         }
390         if (Pstream::master())
391         {
392             outputFilePtr_()<< endl;
393         }
395         if (log_)
396         {
397             Info<< endl;
398         }
399     }
403 // ************************************************************************* //