fieldMinMaxTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fieldMinMax.H"
27 #include "volFields.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const word& fieldName,
35  const word& outputName,
36  const vector& minC,
37  const vector& maxC,
38  const label minProcI,
39  const label maxProcI,
40  const Type& minValue,
41  const Type& maxValue
42 )
43 {
44  OFstream& file = this->file();
45 
46  if (location_)
47  {
48  file<< obr_.time().value();
49 
50  writeTabbed(file, fieldName);
51 
52  file<< token::TAB << minValue
53  << token::TAB << minC;
54 
55  if (Pstream::parRun())
56  {
57  file<< token::TAB << minProcI;
58  }
59 
60  file<< token::TAB << maxValue
61  << token::TAB << maxC;
62 
63  if (Pstream::parRun())
64  {
65  file<< token::TAB << maxProcI;
66  }
67 
68  file<< endl;
69 
70  if (log_) Info
71  << " min(" << outputName << ") = " << minValue
72  << " at location " << minC;
73 
74  if (Pstream::parRun())
75  {
76  if (log_) Info<< " on processor " << minProcI;
77  }
78 
79  if (log_) Info
80  << nl << " max(" << outputName << ") = " << maxValue
81  << " at location " << maxC;
82 
83  if (Pstream::parRun())
84  {
85  if (log_) Info<< " on processor " << maxProcI;
86  }
87  }
88  else
89  {
90  file<< token::TAB << minValue << token::TAB << maxValue;
91 
92  if (log_) Info
93  << " min/max(" << outputName << ") = "
94  << minValue << ' ' << maxValue;
95  }
96 
97  if (log_) Info<< endl;
98 }
99 
100 
101 template<class Type>
103 (
104  const word& fieldName,
105  const modeType& mode
106 )
107 {
109 
110  if (obr_.foundObject<fieldType>(fieldName))
111  {
112  const label procI = Pstream::myProcNo();
113 
114  const fieldType& field = obr_.lookupObject<fieldType>(fieldName);
115  const fvMesh& mesh = field.mesh();
116 
117  const volVectorField::GeometricBoundaryField& CfBoundary =
118  mesh.C().boundaryField();
119 
120  switch (mode)
121  {
122  case mdMag:
123  {
124  const volScalarField magField(mag(field));
125  const volScalarField::GeometricBoundaryField& magFieldBoundary =
126  magField.boundaryField();
127 
128  scalarList minVs(Pstream::nProcs());
129  List<vector> minCs(Pstream::nProcs());
130  label minProcI = findMin(magField);
131  minVs[procI] = magField[minProcI];
132  minCs[procI] = field.mesh().C()[minProcI];
133 
134 
135  labelList maxIs(Pstream::nProcs());
136  scalarList maxVs(Pstream::nProcs());
137  List<vector> maxCs(Pstream::nProcs());
138  label maxProcI = findMax(magField);
139  maxVs[procI] = magField[maxProcI];
140  maxCs[procI] = field.mesh().C()[maxProcI];
141 
142  forAll(magFieldBoundary, patchI)
143  {
144  const scalarField& mfp = magFieldBoundary[patchI];
145  if (mfp.size())
146  {
147  const vectorField& Cfp = CfBoundary[patchI];
148 
149  label minPI = findMin(mfp);
150  if (mfp[minPI] < minVs[procI])
151  {
152  minVs[procI] = mfp[minPI];
153  minCs[procI] = Cfp[minPI];
154  }
155 
156  label maxPI = findMax(mfp);
157  if (mfp[maxPI] > maxVs[procI])
158  {
159  maxVs[procI] = mfp[maxPI];
160  maxCs[procI] = Cfp[maxPI];
161  }
162  }
163  }
164 
165  Pstream::gatherList(minVs);
166  Pstream::gatherList(minCs);
167 
168  Pstream::gatherList(maxVs);
169  Pstream::gatherList(maxCs);
170 
171  if (Pstream::master())
172  {
173  label minI = findMin(minVs);
174  scalar minValue = minVs[minI];
175  const vector& minC = minCs[minI];
176 
177  label maxI = findMax(maxVs);
178  scalar maxValue = maxVs[maxI];
179  const vector& maxC = maxCs[maxI];
180 
181  output
182  (
183  fieldName,
184  word("mag(" + fieldName + ")"),
185  minC,
186  maxC,
187  minI,
188  maxI,
189  minValue,
190  maxValue
191  );
192  }
193  break;
194  }
195  case mdCmpt:
196  {
197  const typename fieldType::GeometricBoundaryField&
198  fieldBoundary = field.boundaryField();
199 
200  List<Type> minVs(Pstream::nProcs());
201  List<vector> minCs(Pstream::nProcs());
202  label minProcI = findMin(field);
203  minVs[procI] = field[minProcI];
204  minCs[procI] = field.mesh().C()[minProcI];
205 
206  Pstream::gatherList(minVs);
207  Pstream::gatherList(minCs);
208 
209  List<Type> maxVs(Pstream::nProcs());
210  List<vector> maxCs(Pstream::nProcs());
211  label maxProcI = findMax(field);
212  maxVs[procI] = field[maxProcI];
213  maxCs[procI] = field.mesh().C()[maxProcI];
214 
215  forAll(fieldBoundary, patchI)
216  {
217  const Field<Type>& fp = fieldBoundary[patchI];
218  if (fp.size())
219  {
220  const vectorField& Cfp = CfBoundary[patchI];
221 
222  label minPI = findMin(fp);
223  if (fp[minPI] < minVs[procI])
224  {
225  minVs[procI] = fp[minPI];
226  minCs[procI] = Cfp[minPI];
227  }
228 
229  label maxPI = findMax(fp);
230  if (fp[maxPI] > maxVs[procI])
231  {
232  maxVs[procI] = fp[maxPI];
233  maxCs[procI] = Cfp[maxPI];
234  }
235  }
236  }
237 
238  Pstream::gatherList(minVs);
239  Pstream::gatherList(minCs);
240 
241  Pstream::gatherList(maxVs);
242  Pstream::gatherList(maxCs);
243 
244  if (Pstream::master())
245  {
246  label minI = findMin(minVs);
247  Type minValue = minVs[minI];
248  const vector& minC = minCs[minI];
249 
250  label maxI = findMax(maxVs);
251  Type maxValue = maxVs[maxI];
252  const vector& maxC = maxCs[maxI];
253 
254  output
255  (
256  fieldName,
257  fieldName,
258  minC,
259  maxC,
260  minI,
261  maxI,
262  minValue,
263  maxValue
264  );
265  }
266  break;
267  }
268  default:
269  {
271  (
272  "Foam::fieldMinMax::calcMinMaxFields"
273  "("
274  "const word&, "
275  "const modeType&"
276  ")"
277  ) << "Unknown min/max mode: " << modeTypeNames_[mode_]
278  << exit(FatalError);
279  }
280  }
281  }
282 }
283 
284 
285 // ************************************************************************* //
Output to file stream.
Definition: OFstream.H:81
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
dimensioned< scalar > mag(const dimensioned< Type > &)
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
label findMax(const ListType &, const label start=0)
Find index of max element (and larger than given element).
A class for handling words, derived from string.
Definition: word.H:59
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
scalar maxValue
messageStream Info
dynamicFvMesh & mesh
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
#define forAll(list, i)
Definition: UList.H:421
label findMin(const ListType &, const label start=0)
Find index of min element (and less than given element).
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
void output(const word &fieldName, const word &outputName, const vector &minC, const vector &maxC, const label minProcI, const label maxProcI, const Type &minValue, const Type &maxValue)
Helper function to write the output.
Generic GeometricField class.
void calcMinMaxFields(const word &fieldName, const modeType &mode)
Calculate the field min/max.
error FatalError