fieldValueDelta.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) 2012-2016 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 "fieldValueDelta.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace functionObjects
34 {
35 namespace fieldValues
36 {
39 }
40 }
41 }
42 
43 template<>
44 const char* Foam::NamedEnum
45 <
47  5
48 >::names[] =
49 {
50  "add",
51  "subtract",
52  "min",
53  "max",
54  "average"
55 };
56 
57 const Foam::NamedEnum
58 <
60  5
62 
63 
64 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
65 
67 (
68  const label i
69 )
70 {
71  const wordList& fields1 = region1Ptr_->fields();
72  const wordList& fields2 = region2Ptr_->fields();
73 
74  DynamicList<word> commonFields(fields1.size());
75  forAll(fields1, i)
76  {
77  label index = findIndex(fields2, fields1[i]);
78  if (index != -1)
79  {
80  commonFields.append(fields1[i]);
81  }
82  }
83 
84  Ostream& os = file();
85 
86  writeHeaderValue(os, "Source1", region1Ptr_->name());
87  writeHeaderValue(os, "Source2", region2Ptr_->name());
88  writeHeaderValue(os, "Operation", operationTypeNames_[operation_]);
89  writeCommented(os, "Time");
90 
91  forAll(commonFields, i)
92  {
93  os << tab << commonFields[i];
94  }
95 
96  os << endl;
97 }
98 
99 
100 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
101 
103 (
104  const word& name,
105  const Time& runTime,
106  const dictionary& dict
107 )
108 :
109  writeFiles(name, runTime, dict, name),
110  operation_(opSubtract),
111  region1Ptr_(NULL),
112  region2Ptr_(NULL)
113 {
114  if (!isA<fvMesh>(obr_))
115  {
117  << "objectRegistry is not an fvMesh" << exit(FatalError);
118  }
119 
120  read(dict);
121  resetName(typeName);
122 }
123 
124 
125 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
126 
128 {}
129 
130 
131 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
132 
134 (
135  const dictionary& dict
136 )
137 {
138  writeFiles::read(dict);
139 
140  region1Ptr_.reset
141  (
143  (
144  name() + ".region1",
145  obr_,
146  dict.subDict("region1"),
147  false
148  ).ptr()
149  );
150  region2Ptr_.reset
151  (
153  (
154  name() + ".region2",
155  obr_,
156  dict.subDict("region2"),
157  false
158  ).ptr()
159  );
160 
161  operation_ = operationTypeNames_.read(dict.lookup("operation"));
162 
163  return true;
164 }
165 
166 
168 {
170 
171  region1Ptr_->write();
172  region2Ptr_->write();
173 
174  if (Pstream::master())
175  {
176  writeTime(file());
177  }
178 
179  Log << type() << " " << name() << " write:" << endl;
180 
181  bool found = false;
182  processFields<scalar>(found);
183  processFields<vector>(found);
184  processFields<sphericalTensor>(found);
185  processFields<symmTensor>(found);
186  processFields<tensor>(found);
187 
188  if (Pstream::master())
189  {
190  file()<< endl;
191  }
192 
193  if (!found)
194  {
195  Log << " none" << endl;
196  }
197  else
198  {
199  Log << endl;
200  }
201 
202  return true;
203 }
204 
205 
207 {
208  return true;
209 }
210 
211 
212 // ************************************************************************* //
virtual void writeFileHeader(const label i)
Output file header information.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const char tab
Definition: Ostream.H:261
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
addToRunTimeSelectionTable(functionObject, fieldValueDelta, dictionary)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual bool read(const dictionary &)
Read from dictionary.
functionObject base class for writing files
Definition: writeFiles.H:56
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:411
static const NamedEnum< operationType, 5 > operationTypeNames_
Operation type names.
Abstract base-class for Time/database function objects.
This function object provides a differencing option between two &#39;field value&#39; function objects...
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Macros for easy insertion into run-time selection tables.
virtual bool read(const dictionary &)
Read optional controls.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
A class for handling words, derived from string.
Definition: word.H:59
defineTypeNameAndDebug(fieldValueDelta, 0)
static autoPtr< fieldValue > New(const word &name, const objectRegistry &obr, const dictionary &dict, const bool output=true)
Return a reference to the selected fieldValue.
Definition: fieldValueNew.C:32
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:297
fieldValueDelta(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
virtual bool write()
Write function.
Definition: writeFiles.C:197
#define Log
Report write to Foam::Info if the local log switch is true.
virtual bool write()
Calculate and write.
bool found
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451