externalCoupledTemplates.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2017 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify i
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "externalCoupled.H"
29 //#include "fvMesh.H"
30 #include "OSspecific.H"
31 #include "Fstream.H"
32 #include "volFields.H"
34 #include "mixedFvPatchFields.H"
37 #include "StringStream.H"
38 #include "globalIndex.H"
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 template<class Type>
43 bool Foam::functionObjects::externalCoupled::readData
44 (
45  const UPtrList<const fvMesh>& meshes,
46  const wordRe& groupName,
47  const word& fieldName
48 )
49 {
50  typedef GeometricField<Type, fvPatchField, volMesh> volFieldType;
51  typedef externalCoupledMixedFvPatchField<Type> patchFieldType;
52 
53  wordList regionNames(meshes.size());
54  forAll(meshes, i)
55  {
56  regionNames[i] = meshes[i].dbDir();
57  }
58 
59  // File only opened on master; contains data for all processors, for all
60  // patchIDs.
61  autoPtr<IFstream> masterFilePtr;
62  if (Pstream::master())
63  {
64  const fileName transferFile
65  (
66  groupDir(commDirectory(), compositeName(regionNames), groupName)
67  / fieldName + ".in"
68  );
69 
70  Log << type() << ": reading data from " << transferFile << endl;
71 
72  masterFilePtr.reset(new IFstream(transferFile));
73 
74  if (!masterFilePtr().good())
75  {
76  FatalIOErrorInFunction(masterFilePtr())
77  << "Cannot open file for region " << compositeName(regionNames)
78  << ", field " << fieldName
79  << exit(FatalIOError);
80  }
81  }
82 
83 
84  label nFound = 0;
85  for (const fvMesh& mesh : meshes)
86  {
87  const volFieldType* vfptr = mesh.findObject<volFieldType>(fieldName);
88 
89  if (!vfptr)
90  {
91  continue;
92  }
93  nFound++;
94 
95  typename volFieldType::Boundary& bf =
96  const_cast<volFieldType*>(vfptr)->boundaryFieldRef();
97 
98  // Get the patches
99  const labelList patchIDs
100  (
101  mesh.boundaryMesh().patchSet
102  (
103  List<wordRe>{groupName}
104  ).sortedToc()
105  );
106 
107  // Handle column-wise reading of patch data. Supports most easy types
108  for (const label patchi : patchIDs)
109  {
110  if (isA<patchFieldType>(bf[patchi]))
111  {
112  // Explicit handling of externalCoupledMixed bcs - they
113  // have specialised reading routines.
114 
115  patchFieldType& pf = refCast<patchFieldType>
116  (
117  bf[patchi]
118  );
119 
120  // Read from master into local stream
121  OStringStream os;
122  readLines
123  (
124  bf[patchi].size(), // number of lines to read
125  masterFilePtr,
126  os
127  );
128 
129  // Pass responsibility for all reading over to bc
130  pf.readData(IStringStream(os.str())());
131 
132  // Update the value from the read coefficient. Bypass any
133  // additional processing by derived type.
134  pf.patchFieldType::evaluate();
135  }
136  else if (isA<mixedFvPatchField<Type>>(bf[patchi]))
137  {
138  mixedFvPatchField<Type>& pf = refCast<mixedFvPatchField<Type>>
139  (
140  bf[patchi]
141  );
142 
143  // Read columns from file for
144  // value, snGrad, refValue, refGrad, valueFraction
145  List<scalarField> data;
146  readColumns
147  (
148  bf[patchi].size(), // number of lines to read
149  4*pTraits<Type>::nComponents+1, // nColumns: 4*Type+1*scalar
150  masterFilePtr,
151  data
152  );
153 
154  // Transfer read data to bc.
155  // Skip value, snGrad
156  direction columni = 2*pTraits<Type>::nComponents;
157 
158  Field<Type>& refValue = pf.refValue();
159  for
160  (
161  direction cmpt = 0;
162  cmpt < pTraits<Type>::nComponents;
163  cmpt++
164  )
165  {
166  refValue.replace(cmpt, data[columni++]);
167  }
168  Field<Type>& refGrad = pf.refGrad();
169  for
170  (
171  direction cmpt = 0;
172  cmpt < pTraits<Type>::nComponents;
173  cmpt++
174  )
175  {
176  refGrad.replace(cmpt, data[columni++]);
177  }
178  pf.valueFraction() = data[columni];
179 
180  // Update the value from the read coefficient. Bypass any
181  // additional processing by derived type.
182  pf.mixedFvPatchField<Type>::evaluate();
183  }
184  else if (isA<fixedGradientFvPatchField<Type>>(bf[patchi]))
185  {
186  fixedGradientFvPatchField<Type>& pf =
187  refCast<fixedGradientFvPatchField<Type>>(bf[patchi]);
188 
189  // Read columns for value and gradient
190  List<scalarField> data;
191  readColumns
192  (
193  bf[patchi].size(), // number of lines to read
194  2*pTraits<Type>::nComponents, // nColumns: Type
195  masterFilePtr,
196  data
197  );
198 
199  // Transfer gradient to bc
200  Field<Type>& gradient = pf.gradient();
201  for
202  (
203  direction cmpt = 0;
204  cmpt < pTraits<Type>::nComponents;
205  cmpt++
206  )
207  {
208  gradient.replace
209  (
210  cmpt,
211  data[pTraits<Type>::nComponents+cmpt]
212  );
213  }
214 
215  // Update the value from the read coefficient. Bypass any
216  // additional processing by derived type.
217  pf.fixedGradientFvPatchField<Type>::evaluate();
218  }
219  else if (isA<fixedValueFvPatchField<Type>>(bf[patchi]))
220  {
221  fixedValueFvPatchField<Type>& pf =
222  refCast<fixedValueFvPatchField<Type>>(bf[patchi]);
223 
224  // Read columns for value only
225  List<scalarField> data;
226  readColumns
227  (
228  bf[patchi].size(), // number of lines to read
229  pTraits<Type>::nComponents, // number of columns to read
230  masterFilePtr,
231  data
232  );
233 
234  // Transfer read value to bc
235  Field<Type> value(bf[patchi].size());
236  for
237  (
238  direction cmpt = 0;
239  cmpt < pTraits<Type>::nComponents;
240  cmpt++
241  )
242  {
243  value.replace(cmpt, data[cmpt]);
244  }
245 
246  pf == value;
247 
248  // Update the value from the read coefficient. Bypass any
249  // additional processing by derived type.
250  pf.fixedValueFvPatchField<Type>::evaluate();
251  }
252  else
253  {
255  << "Unsupported boundary condition " << bf[patchi].type()
256  << " for patch " << bf[patchi].patch().name()
257  << " in region " << mesh.name()
258  << exit(FatalError);
259  }
260 
261  initialisedCoupling_ = true;
262  }
263  }
264 
265  return nFound > 0;
266 }
267 
268 
269 template<class Type>
271 Foam::functionObjects::externalCoupled::gatherAndCombine
272 (
273  const Field<Type>& fld
274 )
275 {
276  // Collect values from all processors
277  List<Field<Type>> gatheredValues(Pstream::nProcs());
278  gatheredValues[Pstream::myProcNo()] = fld;
279  Pstream::gatherList(gatheredValues);
280 
281 
282  auto tresult = tmp<Field<Type>>::New();
283  auto& result = tresult.ref();
284 
285  if (Pstream::master())
286  {
287  // Combine values into single field
288  label globalElemi = 0;
289 
290  forAll(gatheredValues, lsti)
291  {
292  globalElemi += gatheredValues[lsti].size();
293  }
294 
295  result.setSize(globalElemi);
296 
297  globalElemi = 0;
298 
299  forAll(gatheredValues, lsti)
300  {
301  const Field<Type>& sub = gatheredValues[lsti];
302 
303  forAll(sub, elemi)
304  {
305  result[globalElemi++] = sub[elemi];
306  }
307  }
308  }
309 
310  return tresult;
311 }
312 
313 
314 template<class Type>
315 bool Foam::functionObjects::externalCoupled::writeData
316 (
318  const wordRe& groupName,
319  const word& fieldName
320 ) const
321 {
322  typedef GeometricField<Type, fvPatchField, volMesh> volFieldType;
323  typedef externalCoupledMixedFvPatchField<Type> patchFieldType;
324 
325  wordList regionNames(meshes.size());
326  forAll(meshes, i)
327  {
328  regionNames[i] = meshes[i].dbDir();
329  }
330 
331  // File only opened on master; contains data for all processors, for all
332  // patchIDs
333  autoPtr<OFstream> masterFilePtr;
334  if (Pstream::master())
335  {
336  const fileName transferFile
337  (
338  groupDir(commDirectory(), compositeName(regionNames), groupName)
339  / fieldName + ".out"
340  );
341 
342  Log << type() << ": writing data to " << transferFile << endl;
343 
344  masterFilePtr.reset(new OFstream(transferFile));
345 
346  if (!masterFilePtr().good())
347  {
348  FatalIOErrorInFunction(masterFilePtr())
349  << "Cannot open file for region " << compositeName(regionNames)
350  << ", field " << fieldName
351  << exit(FatalIOError);
352  }
353  }
354 
355 
356  bool headerDone = false;
357 
358  label nFound = 0;
359 
360  for (const fvMesh& mesh : meshes)
361  {
362  const volFieldType* vfptr = mesh.findObject<volFieldType>(fieldName);
363 
364  if (!vfptr)
365  {
366  continue;
367  }
368  nFound++;
369 
370  const typename volFieldType::Boundary& bf =
371  vfptr->boundaryField();
372 
373  // Get the patches
374  const labelList patchIDs
375  (
376  mesh.boundaryMesh().patchSet
377  (
378  List<wordRe>{groupName}
379  ).sortedToc()
380  );
381 
382  // Handle column-wise writing of patch data. Supports most easy types
383  for (const label patchi : patchIDs)
384  {
385  const globalIndex globalFaces(bf[patchi].size());
386 
387  if (isA<patchFieldType>(bf[patchi]))
388  {
389  // Explicit handling of externalCoupledMixed bcs - they
390  // have specialised writing routines
391 
392  const patchFieldType& pf = refCast<const patchFieldType>
393  (
394  bf[patchi]
395  );
396  OStringStream os;
397 
398  // Pass responsibility for all writing over to bc
399  pf.writeData(os);
400 
401  // Collect contributions from all processors and output them on
402  // master
403  if (Pstream::master())
404  {
405  // Output master data first
406  if (!headerDone)
407  {
408  pf.writeHeader(masterFilePtr());
409  headerDone = true;
410  }
411  masterFilePtr() << os.str().c_str();
412 
413  for (label proci = 1; proci < Pstream::nProcs(); proci++)
414  {
415  IPstream fromSlave
416  (
418  proci
419  );
420 
421  string str(fromSlave);
422  masterFilePtr() << str.c_str();
423  }
424  }
425  else
426  {
427  OPstream toMaster
428  (
431  );
432 
433  toMaster << os.str();
434  }
435  }
436  else if (isA<mixedFvPatchField<Type>>(bf[patchi]))
437  {
438  const mixedFvPatchField<Type>& pf =
439  refCast<const mixedFvPatchField<Type>>(bf[patchi]);
440 
441  Field<Type> value(gatherAndCombine(pf));
442  Field<Type> snGrad(gatherAndCombine(pf.snGrad()()));
443  Field<Type> refValue(gatherAndCombine(pf.refValue()));
444  Field<Type> refGrad(gatherAndCombine(pf.refGrad()));
445  scalarField valueFraction(gatherAndCombine(pf.valueFraction()));
446 
447  if (Pstream::master())
448  {
449  forAll(refValue, facei)
450  {
451  masterFilePtr()
452  << value[facei] << token::SPACE
453  << snGrad[facei] << token::SPACE
454  << refValue[facei] << token::SPACE
455  << refGrad[facei] << token::SPACE
456  << valueFraction[facei] << nl;
457  }
458  }
459  }
460  else
461  {
462  // Output the value and snGrad
463  Field<Type> value(gatherAndCombine(bf[patchi]));
464  Field<Type> snGrad(gatherAndCombine(bf[patchi].snGrad()()));
465  if (Pstream::master())
466  {
467  forAll(value, facei)
468  {
469  masterFilePtr()
470  << value[facei] << token::SPACE
471  << snGrad[facei] << nl;
472  }
473  }
474  }
475  }
476  }
477 
478  return nFound > 0;
479 }
480 
481 
482 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
volFields.H
Foam::autoPtr::reset
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:158
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:47
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
OSspecific.H
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
Foam::UPstream::masterNo
static constexpr int masterNo() noexcept
Process index of the master.
Definition: UPstream.H:432
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
externalCoupledMixedFvPatchFields.H
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
globalIndex.H
Foam::UPstream::nProcs
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:426
Foam::externalFileCoupler::commDirectory
const fileName & commDirectory() const
Return the file path to the base communications directory.
Definition: externalFileCouplerI.H:42
StringStream.H
Input/output from string buffers.
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::wordRe
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition: wordRe.H:72
externalCoupled.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:59
Foam::externalCoupledMixedFvPatchField
Extends the mixed boundary condition with serialisation functions.
Definition: externalCoupledMixedFvPatchField.H:76
Foam::label
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:62
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::UPtrList
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:63
meshes
PtrList< fvMesh > meshes(regionNames.size())
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::UPstream::commsTypes::scheduled
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
mixedFvPatchFields.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:99
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:444
fixedGradientFvPatchFields.H
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:438
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::isA
const TargetType * isA(const Type &t)
Check if dynamic_cast to TargetType is possible.
Definition: typeInfo.H:197
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:52
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Fstream.H
Input/output from file streams.
Foam::functionObject::type
virtual const word & type() const =0
Runtime type information.
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:102
fixedValueFvPatchFields.H
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::token::SPACE
Space [isspace].
Definition: token.H:112
Foam::direction
uint8_t direction
Definition: direction.H:47
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:375
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::functionObjects::externalCoupled::compositeName
static word compositeName(const wordList &)
Create single name by appending words (in sorted order),.
Definition: externalCoupled.C:351
Log
#define Log
Report write to Foam::Info if the local log switch is true.
Definition: messageStream.H:332