initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / cfdTools / general / fieldSources / timeActivatedExplicitMulticomponentPointSource / timeActivatedExplicitMulticomponentPointSource.C
blob221b7949b4365c9bad6c6753e1cac6b831e579d5
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "timeActivatedExplicitMulticomponentPointSource.H"
28 #include "volFields.H"
30 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
32 Foam::label
33 Foam::timeActivatedExplicitMulticomponentPointSource::carrierFieldId
35     const word& fieldName
38     forAll(carrierFields_, fieldI)
39     {
40         if (carrierFields_[fieldI].name() == fieldName)
41         {
42             return fieldI;
43         }
44     }
46     return -1;
50 void Foam::timeActivatedExplicitMulticomponentPointSource::updateAddressing()
52     forAll(pointSources_, sourceI)
53     {
54         const pointSourceProperties& psp = pointSources_[sourceI];
55         bool foundCell = false;
56         label cid = mesh_.findCell(psp.location());
57         if (cid >= 0)
58         {
59             foundCell = mesh_.pointInCell(psp.location(), cid);
60         }
61         reduce(foundCell, orOp<bool>());
62         if (!foundCell)
63         {
64             label cid = mesh_.findNearestCell(psp.location());
65             if (cid >= 0)
66             {
67                 foundCell = mesh_.pointInCell(psp.location(), cid);
68             }
69         }
70         reduce(foundCell, orOp<bool>());
72         if (!foundCell)
73         {
74             FatalErrorIn
75             (
76                 "timeActivatedExplicitMulticomponentPointSource::"
77                 "updateAddressing()"
78             )   << "Unable to find location " << psp.location() << " in mesh "
79                 << "for source " << psp.name() << nl
80                 << exit(FatalError);
81         }
82         else
83         {
84             cellOwners_[sourceI] = cid;
85         }
87         fieldIds_[sourceI].setSize(psp.fieldData().size());
88         forAll(psp.fieldData(), fieldI)
89         {
90             const word& fieldName = psp.fieldData()[fieldI].first();
91             label cfid = carrierFieldId(fieldName);
92             if (cfid < 0)
93             {
94                 FatalErrorIn
95                 (
96                     "timeActivatedExplicitMulticomponentPointSource::"
97                     "updateAddressing()"
98                 )   << "Unable to find field " << fieldName << " in carrier "
99                     << "fields for source " << psp.name() << nl
100                     << exit(FatalError);
101             }
102             else
103             {
104                 fieldIds_[sourceI][fieldI] = cfid;
105             }
106        }
107     }
111 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
113 Foam::timeActivatedExplicitMulticomponentPointSource::
114 timeActivatedExplicitMulticomponentPointSource
116     const word& name,
117     const fvMesh& mesh,
118     const PtrList<volScalarField>& carrierFields,
119     const dimensionSet& dims
122     IOdictionary
123     (
124         IOobject
125         (
126             name + "Properties",
127             mesh.time().constant(),
128             mesh,
129             IOobject::MUST_READ,
130             IOobject::NO_WRITE
131         )
132     ),
133     name_(name),
134     mesh_(mesh),
135     runTime_(mesh.time()),
136     dimensions_(dims),
137     carrierFields_(carrierFields),
138     active_(lookup("active")),
139     pointSources_(lookup("pointSources")),
140     cellOwners_(pointSources_.size()),
141     fieldIds_(pointSources_.size())
143     // Initialise the field addressing
144     updateAddressing();
148 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
150 Foam::tmp<Foam::DimensionedField<Foam::scalar, Foam::volMesh> >
151 Foam::timeActivatedExplicitMulticomponentPointSource::Su
153     const label fieldI
156     if (mesh_.changing())
157     {
158         updateAddressing();
159     }
161     tmp<DimensionedField<scalar, volMesh> > tSource
162     (
163         new DimensionedField<scalar, volMesh>
164         (
165             IOobject
166             (
167                 name_ + carrierFields_[fieldI].name() + "Su",
168                 runTime_.timeName(),
169                 mesh_,
170                 IOobject::NO_READ,
171                 IOobject::NO_WRITE
172             ),
173             mesh_,
174             dimensionedScalar("zero", dimensions_, 0.0)
175         )
176     );
178     if (active_)
179     {
180         DimensionedField<scalar, volMesh>& sourceField = tSource();
182         const scalarField& V = mesh_.V();
183         const scalar dt = runTime_.deltaT().value();
185         forAll(pointSources_, sourceI)
186         {
187             const pointSourceProperties& psp = pointSources_[sourceI];
189             forAll(fieldIds_[sourceI], i)
190             {
191                 if
192                 (
193                     fieldIds_[sourceI][i] == fieldI
194                 && (runTime_.time().value() >= psp.timeStart())
195                 && (runTime_.time().value() <= psp.timeEnd())
196                 )
197                 {
198                     const label cid = cellOwners_[sourceI];
199                     if (cid >= 0)
200                     {
201                         sourceField[cid] +=
202                             dt*psp.fieldData()[i].second()/V[cid];
203                     }
204                 }
205             }
206         }
207     }
209     return tSource;
213 Foam::tmp<Foam::DimensionedField<Foam::scalar, Foam::volMesh> >
214 Foam::timeActivatedExplicitMulticomponentPointSource::Su()
216     if (mesh_.changing())
217     {
218         updateAddressing();
219     }
221     tmp<DimensionedField<scalar, volMesh> > tSource
222     (
223         new DimensionedField<scalar, volMesh>
224         (
225             IOobject
226             (
227                 name_ + "TotalSu",
228                 runTime_.timeName(),
229                 mesh_,
230                 IOobject::NO_READ,
231                 IOobject::NO_WRITE
232             ),
233             mesh_,
234             dimensionedScalar("zero", dimensions_, 0.0)
235         )
236     );
238     if (active_)
239     {
240         DimensionedField<scalar, volMesh>& sourceField = tSource();
242         const scalarField& V = mesh_.V();
243         const scalar dt = runTime_.deltaT().value();
245         forAll(pointSources_, sourceI)
246         {
247             const pointSourceProperties& psp = pointSources_[sourceI];
249             forAll(fieldIds_[sourceI], i)
250             {
251                 if
252                 (
253                    (runTime_.time().value() >= psp.timeStart())
254                 && (runTime_.time().value() <= psp.timeEnd())
255                 )
256                 {
257                     const label cid = cellOwners_[sourceI];
258                     if (cid >= 0)
259                     {
260                         sourceField[cid] +=
261                             dt*psp.fieldData()[i].second()/V[cid];
262                     }
263                 }
264             }
265         }
266     }
268     return tSource;
272 bool Foam::timeActivatedExplicitMulticomponentPointSource::read()
274     if (regIOobject::read())
275     {
276         lookup("active") >> active_;
277         lookup("pointSources") >> pointSources_;
279         cellOwners_.setSize(pointSources_.size());
281         updateAddressing();
283         return true;
284     }
285     else
286     {
287         return false;
288     }
292 // ************************************************************************* //