Loading...
Searching...
No Matches
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-2024 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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 "SpanStream.H"
37#include "StringStream.H"
38#include "globalIndex.H"
39
40// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41
42template<class Type>
43bool 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
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 (UPstream::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
80 }
81 }
82
83
84 const wordRes patchSelection(Foam::one{}, groupName);
85
86 label nFound = 0;
87
88 for (const fvMesh& mesh : meshes)
89 {
90 auto* vfptr = mesh.getObjectPtr<volFieldType>(fieldName);
91
92 if (!vfptr)
93 {
94 continue;
95 }
96 ++nFound;
97
98 auto& bf = 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 std::string lines;
121 readLines
122 (
123 bf[patchi].size(), // number of lines to read
124 masterFilePtr,
125 lines
126 );
127
128 ISpanStream isstr(lines);
129
130 // Pass responsibility for all reading over to bc
131 pf.readData(isstr);
132
133 // Update the value from the read coefficient. Bypass any
134 // additional processing by derived type.
135 pf.patchFieldType::evaluate();
136 }
137 else if (isA<mixedFvPatchField<Type>>(bf[patchi]))
138 {
139 mixedFvPatchField<Type>& pf = refCast<mixedFvPatchField<Type>>
140 (
141 bf[patchi]
142 );
143
144 // Read columns from file for
145 // value, snGrad, refValue, refGrad, valueFraction
146 List<scalarField> data;
147 readColumns
148 (
149 bf[patchi].size(), // number of lines to read
150 4*pTraits<Type>::nComponents+1, // nColumns: 4*Type+1*scalar
151 masterFilePtr,
152 data
153 );
154
155 // Transfer read data to bc.
156 // Skip value, snGrad
157 direction columni = 2*pTraits<Type>::nComponents;
158
159 Field<Type>& refValue = pf.refValue();
160 for
161 (
162 direction cmpt = 0;
163 cmpt < pTraits<Type>::nComponents;
164 cmpt++
165 )
166 {
167 refValue.replace(cmpt, data[columni++]);
168 }
169 Field<Type>& refGrad = pf.refGrad();
170 for
171 (
172 direction cmpt = 0;
173 cmpt < pTraits<Type>::nComponents;
174 cmpt++
175 )
176 {
177 refGrad.replace(cmpt, data[columni++]);
178 }
179 pf.valueFraction() = data[columni];
180
181 // Update the value from the read coefficient. Bypass any
182 // additional processing by derived type.
183 pf.mixedFvPatchField<Type>::evaluate();
184 }
185 else if (isA<fixedGradientFvPatchField<Type>>(bf[patchi]))
186 {
187 fixedGradientFvPatchField<Type>& pf =
189
190 // Read columns for value and gradient
191 List<scalarField> data;
192 readColumns
193 (
194 bf[patchi].size(), // number of lines to read
195 2*pTraits<Type>::nComponents, // nColumns: Type
196 masterFilePtr,
197 data
198 );
199
200 // Transfer gradient to bc
201 Field<Type>& gradient = pf.gradient();
202 for
203 (
204 direction cmpt = 0;
205 cmpt < pTraits<Type>::nComponents;
206 cmpt++
207 )
208 {
209 gradient.replace
210 (
211 cmpt,
212 data[pTraits<Type>::nComponents+cmpt]
213 );
214 }
215
216 // Update the value from the read coefficient. Bypass any
217 // additional processing by derived type.
218 pf.fixedGradientFvPatchField<Type>::evaluate();
219 }
220 else if (isA<fixedValueFvPatchField<Type>>(bf[patchi]))
221 {
222 fixedValueFvPatchField<Type>& pf =
224
225 // Read columns for value only
226 List<scalarField> data;
227 readColumns
228 (
229 bf[patchi].size(), // number of lines to read
230 pTraits<Type>::nComponents, // number of columns to read
231 masterFilePtr,
232 data
233 );
234
235 // Transfer read value to bc
236 Field<Type> value(bf[patchi].size());
237 for
238 (
239 direction cmpt = 0;
240 cmpt < pTraits<Type>::nComponents;
241 cmpt++
242 )
243 {
244 value.replace(cmpt, data[cmpt]);
245 }
246
247 pf == value;
248
249 // Update the value from the read coefficient. Bypass any
250 // additional processing by derived type.
251 pf.fixedValueFvPatchField<Type>::evaluate();
252 }
253 else
254 {
256 << "Unsupported boundary condition " << bf[patchi].type()
257 << " for patch " << bf[patchi].patch().name()
258 << " in region " << mesh.name()
259 << exit(FatalError);
260 }
261
262 initialisedCoupling_ = true;
263 }
264 }
265
266 return nFound;
267}
268
269
270template<class Type>
271bool Foam::functionObjects::externalCoupled::writeData
272(
273 const UPtrList<const fvMesh>& meshes,
274 const wordRe& groupName,
275 const word& fieldName
276) const
277{
278 typedef GeometricField<Type, fvPatchField, volMesh> volFieldType;
279 typedef externalCoupledMixedFvPatchField<Type> patchFieldType;
280
282 forAll(meshes, i)
283 {
284 regionNames[i] = meshes[i].dbDir();
285 }
286
287 // File only opened on master; contains data for all processors, for all
288 // patchIDs
289 autoPtr<OFstream> masterFilePtr;
290 if (Pstream::master())
291 {
292 const fileName transferFile
293 (
294 groupDir(commDirectory(), compositeName(regionNames), groupName)
295 / fieldName + ".out"
296 );
297
298 Log << type() << ": writing data to " << transferFile << endl;
299
300 masterFilePtr.reset(new OFstream(transferFile));
301
302 if (!masterFilePtr().good())
303 {
304 FatalIOErrorInFunction(masterFilePtr())
305 << "Cannot open file for region " << compositeName(regionNames)
306 << ", field " << fieldName
307 << exit(FatalIOError);
308 }
309 }
310
311
312 const wordRes patchSelection(Foam::one{}, groupName);
313
314 bool headerDone = false;
315 label nFound = 0;
316
317 for (const fvMesh& mesh : meshes)
318 {
319 const auto* vfptr = mesh.getObjectPtr<volFieldType>(fieldName);
320
321 if (!vfptr)
322 {
323 continue;
324 }
325 ++nFound;
326
327 const auto& bf = vfptr->boundaryField();
328
329 // Get the patches
330 const labelList patchIDs
331 (
332 mesh.boundaryMesh().patchSet(patchSelection).sortedToc()
333 );
334
335 // Handle column-wise writing of patch data. Supports most easy types
336 for (const label patchi : patchIDs)
337 {
338 // const globalIndex globalFaces(bf[patchi].size());
339
340 if (isA<patchFieldType>(bf[patchi]))
341 {
342 // Explicit handling of externalCoupledMixed bcs - they
343 // have specialised writing routines
344
345 const patchFieldType& pf = refCast<const patchFieldType>
346 (
347 bf[patchi]
348 );
349 OStringStream os;
350
351 // Pass responsibility for all writing over to bc
352 pf.writeData(os);
353
354 // Collect contributions from all processors and output them on
355 // master
356 if (UPstream::master())
357 {
358 // Output master data first
359 if (!headerDone)
360 {
361 pf.writeHeader(masterFilePtr());
362 headerDone = true;
363 }
364 masterFilePtr() << os.str().c_str();
365
366 for (const int proci : UPstream::subProcs())
367 {
368 string str;
369 IPstream::recv(str, proci);
370
371 masterFilePtr() << str.c_str();
372 }
373 }
374 else
375 {
377 }
378 }
379 else if (isA<mixedFvPatchField<Type>>(bf[patchi]))
380 {
381 const mixedFvPatchField<Type>& pf =
383
384 const globalIndex glob
385 (
386 globalIndex::gatherOnly{},
387 pf.size()
388 );
389
390 Field<Type> value(glob.gather(pf));
391 Field<Type> snGrad(glob.gather(pf.snGrad()()));
392 Field<Type> refValue(glob.gather(pf.refValue()));
393 Field<Type> refGrad(glob.gather(pf.refGrad()));
394 scalarField valueFraction(glob.gather(pf.valueFraction()));
395
396 if (UPstream::master())
397 {
398 forAll(value, facei)
399 {
400 masterFilePtr()
401 << value[facei] << token::SPACE
402 << snGrad[facei] << token::SPACE
403 << refValue[facei] << token::SPACE
404 << refGrad[facei] << token::SPACE
405 << valueFraction[facei] << nl;
406 }
407 }
408 }
409 else
410 {
411 // Output the value and snGrad
412
413 const globalIndex glob
414 (
415 globalIndex::gatherOnly{},
416 bf[patchi].size()
417 );
418
419 Field<Type> value(glob.gather(bf[patchi]));
420 Field<Type> snGrad(glob.gather(bf[patchi].snGrad()()));
421
422 if (UPstream::master())
423 {
424 forAll(value, facei)
425 {
426 masterFilePtr()
427 << value[facei] << token::SPACE
428 << snGrad[facei] << nl;
429 }
430 }
431 }
432 }
433 }
434
435 return nFound;
436}
437
438
439// ************************************************************************* //
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
#define Log
Definition PDRblock.C:28
Input/output streams with (internal or external) character storage.
Input/output from string buffers.
labelList patchIDs
static void recv(Type &value, const int fromProcNo, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm, IOstreamOption::streamFormat fmt=IOstreamOption::BINARY)
Receive and deserialize a value. Uses operator>> for de-serialization.
Definition IPstream.H:80
bool send()
Send buffer contents now and not in destructor [advanced usage]. Returns true on success.
Definition OPstreams.C:84
static constexpr int masterNo() noexcept
Relative rank for the master process - is always 0.
Definition UPstream.H:1691
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static rangeType subProcs(const label communicator=worldComm)
Range of process indices for sub-processes.
Definition UPstream.H:1866
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
const fileName & commDirectory() const
Return the file path to the base communications directory.
virtual const word & type() const =0
Runtime type information.
static word compositeName(const wordList &)
Create single name by appending words (in sorted order), separated by '_'.
@ SPACE
Space [isspace].
Definition token.H:144
dynamicFvMesh & mesh
Foam::PtrList< Foam::fvMesh > meshes(regionNames.size())
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
wordList regionNames
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvcSnGrad.C:40
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< word > wordList
List of word.
Definition fileName.H:60
List< label > labelList
A List of labels.
Definition List.H:62
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
uint8_t direction
Definition direction.H:49
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299