Loading...
Searching...
No Matches
mapDistributeIO.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) 2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
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 "mapDistribute.H"
29#include "dictionary.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
38}
39
40
43 is >> *this;
44}
45
46
47// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
48
50{
52
53 // Treat missing as empty
54 const dictionary* subdict = dict.findDict("transforms");
55
56 if (subdict)
57 {
58 subdict->readIfPresent("elements", transformElements_);
59 subdict->readIfPresent("starts", transformStart_);
60 }
61 else
62 {
63 transformElements_.clear();
64 transformStart_.clear();
65 }
66}
67
68
69void Foam::mapDistribute::writeEntries(Ostream& os) const
70{
72
73 if (transformElements_.size())
74 {
75 os << nl;
76 os.beginBlock("transforms");
77 os.writeEntry("elements", transformElements_);
78 transformStart_.writeEntry("starts", os);
79 os.endBlock();
80 }
81}
82
83
84// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
85
86Foam::Istream& Foam::operator>>(Istream& is, mapDistribute& map)
87{
89
90 is >> static_cast<mapDistributeBase&>(map)
91 >> map.transformElements_ >> map.transformStart_;
92
93 return is;
94}
95
96
97Foam::Ostream& Foam::operator<<(Ostream& os, const mapDistribute& map)
98{
100 << map.transformElements_ << token::NL
101 << map.transformStart_;
102
103 return os;
104}
105
106
107// ************************************************************************* //
bool fatalCheck(const char *operation) const
Check IOstream status for given operation.
Definition IOstream.C:51
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
Class containing processor-to-processor mapping information.
void readDict(const dictionary &dict)
Read entries from dictionary format.
void writeEntries(Ostream &os) const
Write entries in dictionary format.
label comm() const noexcept
The communicator used.
Class containing processor-to-processor mapping information.
void readDict(const dictionary &dict)
Read entries from dictionary format.
mapDistribute() noexcept
Default construct - uses worldComm.
void writeEntries(Ostream &os) const
Write entries in dictionary format.
@ NL
Newline [isspace].
Definition token.H:143
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
Istream & operator>>(Istream &, directionInfo &)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict