1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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
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 * * * * * * * * * * * //
33 Foam::timeActivatedExplicitMulticomponentPointSource::carrierFieldId
38 forAll(carrierFields_, fieldI)
40 if (carrierFields_[fieldI].name() == fieldName)
50 void Foam::timeActivatedExplicitMulticomponentPointSource::updateAddressing()
52 forAll(pointSources_, sourceI)
54 const pointSourceProperties& psp = pointSources_[sourceI];
55 bool foundCell = false;
56 label cid = mesh_.findCell(psp.location());
59 foundCell = mesh_.pointInCell(psp.location(), cid);
61 reduce(foundCell, orOp<bool>());
64 label cid = mesh_.findNearestCell(psp.location());
67 foundCell = mesh_.pointInCell(psp.location(), cid);
70 reduce(foundCell, orOp<bool>());
76 "timeActivatedExplicitMulticomponentPointSource::"
78 ) << "Unable to find location " << psp.location() << " in mesh "
79 << "for source " << psp.name() << nl
84 cellOwners_[sourceI] = cid;
87 fieldIds_[sourceI].setSize(psp.fieldData().size());
88 forAll(psp.fieldData(), fieldI)
90 const word& fieldName = psp.fieldData()[fieldI].first();
91 label cfid = carrierFieldId(fieldName);
96 "timeActivatedExplicitMulticomponentPointSource::"
98 ) << "Unable to find field " << fieldName << " in carrier "
99 << "fields for source " << psp.name() << nl
104 fieldIds_[sourceI][fieldI] = cfid;
111 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
113 Foam::timeActivatedExplicitMulticomponentPointSource::
114 timeActivatedExplicitMulticomponentPointSource
118 const PtrList<volScalarField>& carrierFields,
119 const dimensionSet& dims
127 mesh.time().constant(),
135 runTime_(mesh.time()),
137 carrierFields_(carrierFields),
138 active_(lookup("active")),
139 pointSources_(lookup("pointSources")),
140 cellOwners_(pointSources_.size()),
141 fieldIds_(pointSources_.size())
143 // Initialise the field addressing
148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 Foam::tmp<Foam::DimensionedField<Foam::scalar, Foam::volMesh> >
151 Foam::timeActivatedExplicitMulticomponentPointSource::Su
156 if (mesh_.changing())
161 tmp<DimensionedField<scalar, volMesh> > tSource
163 new DimensionedField<scalar, volMesh>
167 name_ + carrierFields_[fieldI].name() + "Su",
174 dimensionedScalar("zero", dimensions_, 0.0)
180 DimensionedField<scalar, volMesh>& sourceField = tSource();
182 const scalarField& V = mesh_.V();
183 const scalar dt = runTime_.deltaT().value();
185 forAll(pointSources_, sourceI)
187 const pointSourceProperties& psp = pointSources_[sourceI];
189 forAll(fieldIds_[sourceI], i)
193 fieldIds_[sourceI][i] == fieldI
194 && (runTime_.time().value() >= psp.timeStart())
195 && (runTime_.time().value() <= psp.timeEnd())
198 const label cid = cellOwners_[sourceI];
202 dt*psp.fieldData()[i].second()/V[cid];
213 Foam::tmp<Foam::DimensionedField<Foam::scalar, Foam::volMesh> >
214 Foam::timeActivatedExplicitMulticomponentPointSource::Su()
216 if (mesh_.changing())
221 tmp<DimensionedField<scalar, volMesh> > tSource
223 new DimensionedField<scalar, volMesh>
234 dimensionedScalar("zero", dimensions_, 0.0)
240 DimensionedField<scalar, volMesh>& sourceField = tSource();
242 const scalarField& V = mesh_.V();
243 const scalar dt = runTime_.deltaT().value();
245 forAll(pointSources_, sourceI)
247 const pointSourceProperties& psp = pointSources_[sourceI];
249 forAll(fieldIds_[sourceI], i)
253 (runTime_.time().value() >= psp.timeStart())
254 && (runTime_.time().value() <= psp.timeEnd())
257 const label cid = cellOwners_[sourceI];
261 dt*psp.fieldData()[i].second()/V[cid];
272 bool Foam::timeActivatedExplicitMulticomponentPointSource::read()
274 if (regIOobject::read())
276 lookup("active") >> active_;
277 lookup("pointSources") >> pointSources_;
279 cellOwners_.setSize(pointSources_.size());
292 // ************************************************************************* //