ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / lagrangian / intermediate / parcels / Templates / CollidingParcel / CollidingParcelIO.C
blobec94fd37599634e9a63446efbe4a6717e3e77583
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "CollidingParcel.H"
27 #include "IOstreams.H"
28 #include "IOField.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 template<class ParcelType>
33 Foam::string Foam::CollidingParcel<ParcelType>::propHeader =
34     ParcelType::propHeader
35   + " collisionRecordsPairAccessed"
36   + " collisionRecordsPairOrigProcOfOther"
37   + " collisionRecordsPairOrigIdOfOther"
38   + " (collisionRecordsPairData)"
39   + " collisionRecordsWallAccessed"
40   + " collisionRecordsWallPRel"
41   + " (collisionRecordsWallData)";
44 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
46 template<class ParcelType>
47 Foam::CollidingParcel<ParcelType>::CollidingParcel
49     const polyMesh& mesh,
50     Istream& is,
51     bool readFields
54     ParcelType(mesh, is, readFields),
55     collisionRecords_()
57     if (readFields)
58     {
59         is >> collisionRecords_;
60     }
62     // Check state of Istream
63     is.check
64     (
65         "CollidingParcel<ParcelType>::Collisions"
66         "(const polyMesh&, Istream&, bool)"
67     );
71 template<class ParcelType>
72 template<class CloudType>
73 void Foam::CollidingParcel<ParcelType>::readFields(CloudType& c)
75     if (!c.size())
76     {
77         return;
78     }
80     ParcelType::readFields(c);
82     labelFieldCompactIOField collisionRecordsPairAccessed
83     (
84         c.fieldIOobject("collisionRecordsPairAccessed", IOobject::MUST_READ)
85     );
86     c.checkFieldFieldIOobject(c, collisionRecordsPairAccessed);
88     labelFieldCompactIOField collisionRecordsPairOrigProcOfOther
89     (
90         c.fieldIOobject
91         (
92             "collisionRecordsPairOrigProcOfOther",
93             IOobject::MUST_READ
94         )
95     );
96     c.checkFieldFieldIOobject(c, collisionRecordsPairOrigProcOfOther);
98     labelFieldCompactIOField collisionRecordsPairOrigIdOfOther
99     (
100         c.fieldIOobject
101         (
102             "collisionRecordsPairOrigIdOfOther",
103             IOobject::MUST_READ
104         )
105     );
106     c.checkFieldFieldIOobject(c, collisionRecordsPairOrigProcOfOther);
108     pairDataFieldCompactIOField collisionRecordsPairData
109     (
110         c.fieldIOobject("collisionRecordsPairData", IOobject::MUST_READ)
111     );
112     c.checkFieldFieldIOobject(c, collisionRecordsPairData);
114     labelFieldCompactIOField collisionRecordsWallAccessed
115     (
116         c.fieldIOobject("collisionRecordsWallAccessed", IOobject::MUST_READ)
117     );
118     c.checkFieldFieldIOobject(c, collisionRecordsWallAccessed);
120     vectorFieldCompactIOField collisionRecordsWallPRel
121     (
122         c.fieldIOobject("collisionRecordsWallPRel", IOobject::MUST_READ)
123     );
124     c.checkFieldFieldIOobject(c, collisionRecordsWallPRel);
126     wallDataFieldCompactIOField collisionRecordsWallData
127     (
128         c.fieldIOobject("collisionRecordsWallData", IOobject::MUST_READ)
129     );
130     c.checkFieldFieldIOobject(c, collisionRecordsWallData);
132     label i = 0;
134     forAllIter(typename CloudType, c, iter)
135     {
136         CollidingParcel<ParcelType>& p = iter();
138         p.collisionRecords_ = collisionRecordList
139         (
140             collisionRecordsPairAccessed[i],
141             collisionRecordsPairOrigProcOfOther[i],
142             collisionRecordsPairOrigIdOfOther[i],
143             collisionRecordsPairData[i],
144             collisionRecordsWallAccessed[i],
145             collisionRecordsWallPRel[i],
146             collisionRecordsWallData[i]
147         );
149         i++;
150     }
154 template<class ParcelType>
155 template<class CloudType>
156 void Foam::CollidingParcel<ParcelType>::writeFields(const CloudType& c)
158     ParcelType::writeFields(c);
160     label np =  c.size();
162     labelFieldCompactIOField collisionRecordsPairAccessed
163     (
164         c.fieldIOobject("collisionRecordsPairAccessed", IOobject::NO_READ),
165         np
166     );
167     labelFieldCompactIOField collisionRecordsPairOrigProcOfOther
168     (
169         c.fieldIOobject
170         (
171             "collisionRecordsPairOrigProcOfOther",
172             IOobject::NO_READ
173         ),
174         np
175     );
176     labelFieldCompactIOField collisionRecordsPairOrigIdOfOther
177     (
178         c.fieldIOobject("collisionRecordsPairOrigIdOfOther", IOobject::NO_READ),
179         np
180     );
181     pairDataFieldCompactIOField collisionRecordsPairData
182     (
183         c.fieldIOobject("collisionRecordsPairData", IOobject::NO_READ),
184         np
185     );
186     labelFieldCompactIOField collisionRecordsWallAccessed
187     (
188         c.fieldIOobject("collisionRecordsWallAccessed", IOobject::NO_READ),
189         np
190     );
191     vectorFieldCompactIOField collisionRecordsWallPRel
192     (
193         c.fieldIOobject("collisionRecordsWallPRel", IOobject::NO_READ),
194         np
195     );
196     wallDataFieldCompactIOField collisionRecordsWallData
197     (
198         c.fieldIOobject("collisionRecordsWallData", IOobject::NO_READ),
199         np
200     );
202     label i = 0;
204     forAllConstIter(typename CloudType, c, iter)
205     {
206         const CollidingParcel<ParcelType>& p = iter();
208         collisionRecordsPairAccessed[i] = p.collisionRecords().pairAccessed();
209         collisionRecordsPairOrigProcOfOther[i] =
210             p.collisionRecords().pairOrigProcOfOther();
211         collisionRecordsPairOrigIdOfOther[i] =
212             p.collisionRecords().pairOrigIdOfOther();
213         collisionRecordsPairData[i] = p.collisionRecords().pairData();
214         collisionRecordsWallAccessed[i] = p.collisionRecords().wallAccessed();
215         collisionRecordsWallPRel[i] = p.collisionRecords().wallPRel();
216         collisionRecordsWallData[i] = p.collisionRecords().wallData();
218         i++;
219     }
221     collisionRecordsPairAccessed.write();
222     collisionRecordsPairOrigProcOfOther.write();
223     collisionRecordsPairOrigIdOfOther.write();
224     collisionRecordsPairData.write();
225     collisionRecordsWallAccessed.write();
226     collisionRecordsWallPRel.write();
227     collisionRecordsWallData.write();
231 // * * * * * * * * * * * * * * * IOstream Operators  * * * * * * * * * * * * //
233 template<class ParcelType>
234 Foam::Ostream& Foam::operator<<
236     Ostream& os,
237     const CollidingParcel<ParcelType>& p
240     if (os.format() == IOstream::ASCII)
241     {
242         os  << static_cast<const ParcelType&>(p)
243             << token::SPACE << p.collisionRecords();
244     }
245     else
246     {
247         os  << static_cast<const ParcelType&>(p)
248             << p.collisionRecords();
249     }
251     // Check state of Ostream
252     os.check
253     (
254         "Ostream& operator<<(Ostream&, const CollidingParcel<ParcelType>&)"
255     );
257     return os;
261 // ************************************************************************* //