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-2021 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 "OSspecific.H"
30 #include "Fstream.H"
31 #include "volFields.H"
33 #include "mixedFvPatchFields.H"
36 #include "StringStream.H"
37 #include "globalIndex.H"
38 
39 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40 
41 template<class Type>
42 bool Foam::functionObjects::externalCoupled::readData
43 (
44  const UPtrList<const fvMesh>& meshes,
45  const wordRe& groupName,
46  const word& fieldName
47 )
48 {
49  typedef GeometricField<Type, fvPatchField, volMesh> volFieldType;
50  typedef externalCoupledMixedFvPatchField<Type> patchFieldType;
51 
52  wordList regionNames(meshes.size());
53  forAll(meshes, i)
54  {
55  regionNames[i] = meshes[i].dbDir();
56  }
57 
58  // File only opened on master; contains data for all processors, for all
59  // patchIDs.
60  autoPtr<IFstream> masterFilePtr;
61  if (Pstream::master())
62  {
63  const fileName transferFile
64  (
65  groupDir(commDirectory(), compositeName(regionNames), groupName)
66  / fieldName + ".in"
67  );
68 
69  Log << type() << ": reading data from " << transferFile << endl;
70 
71  masterFilePtr.reset(new IFstream(transferFile));
72 
73  if (!masterFilePtr().good())
74  {
75  FatalIOErrorInFunction(masterFilePtr())
76  << "Cannot open file for region " << compositeName(regionNames)
77  << ", field " << fieldName
78  << exit(FatalIOError);
79  }
80  }
81 
82 
83  const wordRes patchSelection(one{}, groupName);
84 
85  label nFound = 0;
86 
87  for (const fvMesh& mesh : meshes)
88  {
89  const volFieldType* vfptr = mesh.findObject<volFieldType>(fieldName);
90 
91  if (!vfptr)
92  {
93  continue;
94  }
95  ++nFound;
96 
97  typename volFieldType::Boundary& bf =
98  const_cast<volFieldType*>(vfptr)->boundaryFieldRef();
99 
100  // Get the patches
101  const labelList patchIDs
102  (
103  mesh.boundaryMesh().patchSet(patchSelection).sortedToc()
104  );
105 
106  // Handle column-wise reading of patch data. Supports most easy types
107  for (const label patchi : patchIDs)
108  {
109  if (isA<patchFieldType>(bf[patchi]))
110  {
111  // Explicit handling of externalCoupledMixed bcs - they
112  // have specialised reading routines.
113 
114  patchFieldType& pf = refCast<patchFieldType>
115  (
116  bf[patchi]
117  );
118 
119  // Read from master into local stream
120  OStringStream os;
121  readLines
122  (
123  bf[patchi].size(), // number of lines to read
124  masterFilePtr,
125  os
126  );
127 
128  // Pass responsibility for all reading over to bc
129  pf.readData(IStringStream(os.str())());
130 
131  // Update the value from the read coefficient. Bypass any
132  // additional processing by derived type.
133  pf.patchFieldType::evaluate();
134  }
135  else if (isA<mixedFvPatchField<Type>>(bf[patchi]))
136  {
137  mixedFvPatchField<Type>& pf = refCast<mixedFvPatchField<Type>>
138  (
139  bf[patchi]
140  );
141 
142  // Read columns from file for
143  // value, snGrad, refValue, refGrad, valueFraction
144  List<scalarField> data;
145  readColumns
146  (
147  bf[patchi].size(), // number of lines to read
148  4*pTraits<Type>::nComponents+1, // nColumns: 4*Type+1*scalar
149  masterFilePtr,
150  data
151  );
152 
153  // Transfer read data to bc.
154  // Skip value, snGrad
155  direction columni = 2*pTraits<Type>::nComponents;
156 
157  Field<Type>& refValue = pf.refValue();
158  for
159  (
160  direction cmpt = 0;
161  cmpt < pTraits<Type>::nComponents;
162  cmpt++
163  )
164  {
165  refValue.replace(cmpt, data[columni++]);
166  }
167  Field<Type>& refGrad = pf.refGrad();
168  for
169  (
170  direction cmpt = 0;
171  cmpt < pTraits<Type>::nComponents;
172  cmpt++
173  )
174  {
175  refGrad.replace(cmpt, data[columni++]);
176  }
177  pf.valueFraction() = data[columni];
178 
179  // Update the value from the read coefficient. Bypass any
180  // additional processing by derived type.
181  pf.mixedFvPatchField<Type>::evaluate();
182  }
183  else if (isA<fixedGradientFvPatchField<Type>>(bf[patchi]))
184  {
185  fixedGradientFvPatchField<Type>& pf =
186  refCast<fixedGradientFvPatchField<Type>>(bf[patchi]);
187 
188  // Read columns for value and gradient
189  List<scalarField> data;
190  readColumns
191  (
192  bf[patchi].size(), // number of lines to read
193  2*pTraits<Type>::nComponents, // nColumns: Type
194  masterFilePtr,
195  data
196  );
197 
198  // Transfer gradient to bc
199  Field<Type>& gradient = pf.gradient();
200  for
201  (
202  direction cmpt = 0;
203  cmpt < pTraits<Type>::nComponents;
204  cmpt++
205  )
206  {
207  gradient.replace
208  (
209  cmpt,
210  data[pTraits<Type>::nComponents+cmpt]
211  );
212  }
213 
214  // Update the value from the read coefficient. Bypass any
215  // additional processing by derived type.
216  pf.fixedGradientFvPatchField<Type>::evaluate();
217  }
218  else if (isA<fixedValueFvPatchField<Type>>(bf[patchi]))
219  {
220  fixedValueFvPatchField<Type>& pf =
221  refCast<fixedValueFvPatchField<Type>>(bf[patchi]);
222 
223  // Read columns for value only
224  List<scalarField> data;
225  readColumns
226  (
227  bf[patchi].size(), // number of lines to read
228  pTraits<Type>::nComponents, // number of columns to read
229  masterFilePtr,
230  data
231  );
232 
233  // Transfer read value to bc
234  Field<Type> value(bf[patchi].size());
235  for
236  (
237  direction cmpt = 0;
238  cmpt < pTraits<Type>::nComponents;
239  cmpt++
240  )
241  {
242  value.replace(cmpt, data[cmpt]);
243  }
244 
245  pf == value;
246 
247  // Update the value from the read coefficient. Bypass any
248  // additional processing by derived type.
249  pf.fixedValueFvPatchField<Type>::evaluate();
250  }
251  else
252  {
254  << "Unsupported boundary condition " << bf[patchi].type()
255  << " for patch " << bf[patchi].patch().name()
256  << " in region " << mesh.name()
257  << exit(FatalError);
258  }
259 
260  initialisedCoupling_ = true;
261  }
262  }
263 
264  return nFound;
265 }
266 
267 
268 template<class Type>
270 Foam::functionObjects::externalCoupled::gatherAndCombine
271 (
272  const Field<Type>& fld
273 )
274 {
275  // Collect values from all processors
276  List<Field<Type>> gatheredValues(Pstream::nProcs());
277  gatheredValues[Pstream::myProcNo()] = fld;
278  Pstream::gatherList(gatheredValues);
279 
280 
281  auto tresult = tmp<Field<Type>>::New();
282  auto& result = tresult.ref();
283 
284  if (Pstream::master())
285  {
286  // Combine values into single field
287  label globalElemi = 0;
288 
289  forAll(gatheredValues, lsti)
290  {
291  globalElemi += gatheredValues[lsti].size();
292  }
293 
294  result.setSize(globalElemi);
295 
296  globalElemi = 0;
297 
298  forAll(gatheredValues, lsti)
299  {
300  const Field<Type>& sub = gatheredValues[lsti];
301 
302  forAll(sub, elemi)
303  {
304  result[globalElemi++] = sub[elemi];
305  }
306  }
307  }
308 
309  return tresult;
310 }
311 
312 
313 template<class Type>
314 bool Foam::functionObjects::externalCoupled::writeData
315 (
317  const wordRe& groupName,
318  const word& fieldName
319 ) const
320 {
321  typedef GeometricField<Type, fvPatchField, volMesh> volFieldType;
322  typedef externalCoupledMixedFvPatchField<Type> patchFieldType;
323 
324  wordList regionNames(meshes.size());
325  forAll(meshes, i)
326  {
327  regionNames[i] = meshes[i].dbDir();
328  }
329 
330  // File only opened on master; contains data for all processors, for all
331  // patchIDs
332  autoPtr<OFstream> masterFilePtr;
333  if (Pstream::master())
334  {
335  const fileName transferFile
336  (
337  groupDir(commDirectory(), compositeName(regionNames), groupName)
338  / fieldName + ".out"
339  );
340 
341  Log << type() << ": writing data to " << transferFile << endl;
342 
343  masterFilePtr.reset(new OFstream(transferFile));
344 
345  if (!masterFilePtr().good())
346  {
347  FatalIOErrorInFunction(masterFilePtr())
348  << "Cannot open file for region " << compositeName(regionNames)
349  << ", field " << fieldName
350  << exit(FatalIOError);
351  }
352  }
353 
354 
355  const wordRes patchSelection(one{}, groupName);
356 
357  bool headerDone = false;
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(patchSelection).sortedToc()
377  );
378 
379  // Handle column-wise writing of patch data. Supports most easy types
380  for (const label patchi : patchIDs)
381  {
382  const globalIndex globalFaces(bf[patchi].size());
383 
384  if (isA<patchFieldType>(bf[patchi]))
385  {
386  // Explicit handling of externalCoupledMixed bcs - they
387  // have specialised writing routines
388 
389  const patchFieldType& pf = refCast<const patchFieldType>
390  (
391  bf[patchi]
392  );
393  OStringStream os;
394 
395  // Pass responsibility for all writing over to bc
396  pf.writeData(os);
397 
398  // Collect contributions from all processors and output them on
399  // master
400  if (Pstream::master())
401  {
402  // Output master data first
403  if (!headerDone)
404  {
405  pf.writeHeader(masterFilePtr());
406  headerDone = true;
407  }
408  masterFilePtr() << os.str().c_str();
409 
410  for (const int proci : Pstream::subProcs())
411  {
412  IPstream fromSlave
413  (
415  proci
416  );
417 
418  string str(fromSlave);
419  masterFilePtr() << str.c_str();
420  }
421  }
422  else
423  {
424  OPstream toMaster
425  (
428  );
429 
430  toMaster << os.str();
431  }
432  }
433  else if (isA<mixedFvPatchField<Type>>(bf[patchi]))
434  {
435  const mixedFvPatchField<Type>& pf =
436  refCast<const mixedFvPatchField<Type>>(bf[patchi]);
437 
438  Field<Type> value(gatherAndCombine(pf));
439  Field<Type> snGrad(gatherAndCombine(pf.snGrad()()));
440  Field<Type> refValue(gatherAndCombine(pf.refValue()));
441  Field<Type> refGrad(gatherAndCombine(pf.refGrad()));
442  scalarField valueFraction(gatherAndCombine(pf.valueFraction()));
443 
444  if (Pstream::master())
445  {
446  forAll(refValue, facei)
447  {
448  masterFilePtr()
449  << value[facei] << token::SPACE
450  << snGrad[facei] << token::SPACE
451  << refValue[facei] << token::SPACE
452  << refGrad[facei] << token::SPACE
453  << valueFraction[facei] << nl;
454  }
455  }
456  }
457  else
458  {
459  // Output the value and snGrad
460  Field<Type> value(gatherAndCombine(bf[patchi]));
461  Field<Type> snGrad(gatherAndCombine(bf[patchi].snGrad()()));
462  if (Pstream::master())
463  {
464  forAll(value, facei)
465  {
466  masterFilePtr()
467  << value[facei] << token::SPACE
468  << snGrad[facei] << nl;
469  }
470  }
471  }
472  }
473  }
474 
475  return nFound;
476 }
477 
478 
479 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
volFields.H
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...
Log
#define Log
Definition: PDRblock.C:35
Foam::UPstream::masterNo
static constexpr int masterNo() noexcept
Process index of the master (always 0)
Definition: UPstream.H:451
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
externalCoupledMixedFvPatchFields.H
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
globalIndex.H
Foam::externalFileCoupler::commDirectory
const fileName & commDirectory() const
Return the file path to the base communications directory.
Definition: externalFileCouplerI.H:42
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
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:369
Foam::wordRe
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition: wordRe.H:80
externalCoupled.H
regionNames
wordList regionNames
Definition: getAllRegionOptions.H:37
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:62
Foam::externalCoupledMixedFvPatchField
Extends the mixed boundary condition with serialisation functions.
Definition: externalCoupledMixedFvPatchField.H:76
Foam::UPstream::subProcs
static rangeType subProcs(const label communicator=worldComm)
Range of process indices for sub-processes.
Definition: UPstream.H:515
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
meshes
Foam::PtrList< Foam::fvMesh > meshes(regionNames.size())
Foam::UPtrList
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:62
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::FatalError
error FatalError
os
OBJstream os(runTime.globalPath()/outputName)
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:53
fixedGradientFvPatchFields.H
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::reset
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
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:453
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::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Fstream.H
Foam::UPstream::commsTypes::scheduled
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: BitOps.H:63
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:125
Foam::direction
uint8_t direction
Definition: direction.H:52
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::functionObjects::externalCoupled::compositeName
static word compositeName(const wordList &)
Definition: externalCoupled.C:345
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:445