FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

sampledSetsTemplates.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "sampledSets.H"
00027 #include <finiteVolume/volFields.H>
00028 #include <OpenFOAM/ListListOps.H>
00029 
00030 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00031 
00032 template <class Type>
00033 Foam::sampledSets::volFieldSampler<Type>::volFieldSampler
00034 (
00035     const word& interpolationScheme,
00036     const GeometricField<Type, fvPatchField, volMesh>& field,
00037     const PtrList<sampledSet>& samplers
00038 )
00039 :
00040     List<Field<Type> >(samplers.size()),
00041     name_(field.name())
00042 {
00043     autoPtr<interpolation<Type> > interpolator
00044     (
00045         interpolation<Type>::New(interpolationScheme, field)
00046     );
00047 
00048     forAll(samplers, seti)
00049     {
00050         Field<Type>& values = this->operator[](seti);
00051         const sampledSet& samples = samplers[seti];
00052 
00053         values.setSize(samples.size());
00054         forAll(samples, samplei)
00055         {
00056             const point& samplePt = samples[samplei];
00057             label celli = samples.cells()[samplei];
00058             label facei = samples.faces()[samplei];
00059 
00060             if (celli == -1 && facei == -1)
00061             {
00062                 // Special condition for illegal sampling points
00063                 values[samplei] = pTraits<Type>::max;
00064             }
00065             else
00066             {
00067                 values[samplei] = interpolator().interpolate
00068                 (
00069                     samplePt,
00070                     celli,
00071                     facei
00072                 );
00073             }
00074         }
00075     }
00076 }
00077 
00078 
00079 template <class Type>
00080 Foam::sampledSets::volFieldSampler<Type>::volFieldSampler
00081 (
00082     const GeometricField<Type, fvPatchField, volMesh>& field,
00083     const PtrList<sampledSet>& samplers
00084 )
00085 :
00086     List<Field<Type> >(samplers.size()),
00087     name_(field.name())
00088 {
00089     forAll(samplers, seti)
00090     {
00091         Field<Type>& values = this->operator[](seti);
00092         const sampledSet& samples = samplers[seti];
00093 
00094         values.setSize(samples.size());
00095         forAll(samples, samplei)
00096         {
00097             label celli = samples.cells()[samplei];
00098 
00099             if (celli ==-1)
00100             {
00101                 values[samplei] = pTraits<Type>::max;
00102             }
00103             else
00104             {
00105                 values[samplei] = field[celli];
00106             }
00107         }
00108     }
00109 }
00110 
00111 
00112 template <class Type>
00113 Foam::sampledSets::volFieldSampler<Type>::volFieldSampler
00114 (
00115     const List<Field<Type> >& values,
00116     const word& name
00117 )
00118 :
00119     List<Field<Type> >(values),
00120     name_(name)
00121 {}
00122 
00123 
00124 template<class Type>
00125 Foam::label Foam::sampledSets::grep
00126 (
00127     fieldGroup<Type>& fieldList,
00128     const wordList& fieldTypes
00129 ) const
00130 {
00131     fieldList.setSize(fieldNames_.size());
00132     label nFields = 0;
00133 
00134     forAll(fieldNames_, fieldi)
00135     {
00136         if
00137         (
00138             fieldTypes[fieldi]
00139          == GeometricField<Type, fvPatchField, volMesh>::typeName
00140         )
00141         {
00142             fieldList[nFields] = fieldNames_[fieldi];
00143             nFields++;
00144         }
00145     }
00146 
00147     fieldList.setSize(nFields);
00148 
00149     return nFields;
00150 }
00151 
00152 
00153 template<class Type>
00154 void Foam::sampledSets::writeSampleFile
00155 (
00156     const coordSet& masterSampleSet,
00157     const PtrList<volFieldSampler<Type> >& masterFields,
00158     const label seti,
00159     const fileName& timeDir,
00160     const writer<Type>& formatter
00161 )
00162 {
00163     wordList valueSetNames(masterFields.size());
00164     List<const Field<Type>*> valueSets(masterFields.size());
00165 
00166     forAll(masterFields, fieldi)
00167     {
00168         valueSetNames[fieldi] = masterFields[fieldi].name();
00169         valueSets[fieldi] = &masterFields[fieldi][seti];
00170     }
00171 
00172     fileName fName
00173     (
00174         timeDir/formatter.getFileName(masterSampleSet, valueSetNames)
00175     );
00176 
00177     OFstream ofs(fName);
00178     if (ofs.opened())
00179     {
00180         formatter.write
00181         (
00182             masterSampleSet,
00183             valueSetNames,
00184             valueSets,
00185             ofs
00186         );
00187     }
00188     else
00189     {
00190         WarningIn
00191         (
00192             "void Foam::sampledSets::writeSampleFile"
00193             "("
00194                 "const coordSet&, "
00195                 "const PtrList<volFieldSampler<Type> >&, "
00196                 "const label, "
00197                 "const fileName&, "
00198                 "const writer<Type>&"
00199             ")"
00200         )   << "File " << ofs.name() << " could not be opened. "
00201             << "No data will be written" << endl;
00202     }
00203 }
00204 
00205 
00206 template<class T>
00207 void Foam::sampledSets::combineSampledValues
00208 (
00209     const PtrList<volFieldSampler<T> >& sampledFields,
00210     const labelListList& indexSets,
00211     PtrList<volFieldSampler<T> >& masterFields
00212 )
00213 {
00214     forAll(sampledFields, fieldi)
00215     {
00216         List<Field<T> > masterValues(indexSets.size());
00217 
00218         forAll(indexSets, seti)
00219         {
00220             // Collect data from all processors
00221             List<Field<T> > gatheredData(Pstream::nProcs());
00222             gatheredData[Pstream::myProcNo()] = sampledFields[fieldi][seti];
00223             Pstream::gatherList(gatheredData);
00224 
00225             if (Pstream::master())
00226             {
00227                 Field<T> allData
00228                 (
00229                     ListListOps::combine<Field<T> >
00230                     (
00231                         gatheredData,
00232                         Foam::accessOp<Field<T> >()
00233                     )
00234                 );
00235 
00236                 masterValues[seti] = UIndirectList<T>
00237                 (
00238                     allData,
00239                     indexSets[seti]
00240                 )();
00241             }
00242         }
00243 
00244         masterFields.set
00245         (
00246             fieldi,
00247             new volFieldSampler<T>
00248             (
00249                 masterValues,
00250                 sampledFields[fieldi].name()
00251             )
00252         );
00253     }
00254 }
00255 
00256 
00257 template<class Type>
00258 void Foam::sampledSets::sampleAndWrite
00259 (
00260     fieldGroup<Type>& fields
00261 )
00262 {
00263     if (fields.size())
00264     {
00265         bool interpolate = interpolationScheme_ != "cell";
00266 
00267         // Create or use existing writer
00268         if (fields.formatter.empty())
00269         {
00270             fields.formatter = writer<Type>::New(writeFormat_);
00271         }
00272 
00273         // Storage for interpolated values
00274         PtrList<volFieldSampler<Type> > sampledFields(fields.size());
00275 
00276         forAll(fields, fieldi)
00277         {
00278             if (Pstream::master() && verbose_)
00279             {
00280                 Pout<< "sampledSets::sampleAndWrite: "
00281                     << fields[fieldi] << endl;
00282             }
00283 
00284             if (loadFromFiles_)
00285             {
00286                 GeometricField<Type, fvPatchField, volMesh> vf
00287                 (
00288                     IOobject
00289                     (
00290                         fields[fieldi],
00291                         mesh_.time().timeName(),
00292                         mesh_,
00293                         IOobject::MUST_READ,
00294                         IOobject::NO_WRITE,
00295                         false
00296                     ),
00297                     mesh_
00298                 );
00299 
00300                 if (interpolate)
00301                 {
00302                     sampledFields.set
00303                     (
00304                         fieldi,
00305                         new volFieldSampler<Type>
00306                         (
00307                             interpolationScheme_,
00308                             vf,
00309                             *this
00310                         )
00311                     );
00312                 }
00313                 else
00314                 {
00315                     sampledFields.set
00316                     (
00317                         fieldi,
00318                         new volFieldSampler<Type>(vf, *this)
00319                     );
00320                 }
00321             }
00322             else
00323             {
00324                 if (interpolate)
00325                 {
00326                     sampledFields.set
00327                     (
00328                         fieldi,
00329                         new volFieldSampler<Type>
00330                         (
00331                             interpolationScheme_,
00332                             mesh_.lookupObject
00333                             <GeometricField<Type, fvPatchField, volMesh> >
00334                             (fields[fieldi]),
00335                             *this
00336                         )
00337                     );
00338                 }
00339                 else
00340                 {
00341                     sampledFields.set
00342                     (
00343                         fieldi,
00344                         new volFieldSampler<Type>
00345                         (
00346                             mesh_.lookupObject
00347                             <GeometricField<Type, fvPatchField, volMesh> >
00348                             (fields[fieldi]),
00349                             *this
00350                         )
00351                     );
00352                 }
00353             }
00354         }
00355 
00356         // Combine sampled fields from processors.
00357         // Note: only master results are valid
00358 
00359         PtrList<volFieldSampler<Type> > masterFields(sampledFields.size());
00360         combineSampledValues(sampledFields, indexSets_, masterFields);
00361 
00362         if (Pstream::master())
00363         {
00364             forAll(masterSampledSets_, seti)
00365             {
00366                 writeSampleFile
00367                 (
00368                     masterSampledSets_[seti],
00369                     masterFields,
00370                     seti,
00371                     outputPath_/mesh_.time().timeName(),
00372                     fields.formatter()
00373                 );
00374             }
00375         }
00376     }
00377 }
00378 
00379 
00380 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines