00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include <finiteVolume/surfaceInterpolate.H>
00027
00028
00029
00030 namespace Foam
00031 {
00032
00033
00034
00035 namespace fvc
00036 {
00037
00038
00039
00040
00041 template<class Type>
00042 tmp<surfaceInterpolationScheme<Type> > scheme
00043 (
00044 const surfaceScalarField& faceFlux,
00045 Istream& streamData
00046 )
00047 {
00048 return surfaceInterpolationScheme<Type>::New
00049 (
00050 faceFlux.mesh(),
00051 faceFlux,
00052 streamData
00053 );
00054 }
00055
00056
00057
00058 template<class Type>
00059 tmp<surfaceInterpolationScheme<Type> > scheme
00060 (
00061 const surfaceScalarField& faceFlux,
00062 const word& name
00063 )
00064 {
00065 return surfaceInterpolationScheme<Type>::New
00066 (
00067 faceFlux.mesh(),
00068 faceFlux,
00069 faceFlux.mesh().interpolationScheme(name)
00070 );
00071 }
00072
00073
00074
00075 template<class Type>
00076 tmp<surfaceInterpolationScheme<Type> > scheme
00077 (
00078 const fvMesh& mesh,
00079 Istream& streamData
00080 )
00081 {
00082 return surfaceInterpolationScheme<Type>::New
00083 (
00084 mesh,
00085 streamData
00086 );
00087 }
00088
00089
00090
00091 template<class Type>
00092 tmp<surfaceInterpolationScheme<Type> > scheme
00093 (
00094 const fvMesh& mesh,
00095 const word& name
00096 )
00097 {
00098 return surfaceInterpolationScheme<Type>::New
00099 (
00100 mesh,
00101 mesh.interpolationScheme(name)
00102 );
00103 }
00104
00105
00106
00107 template<class Type>
00108 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
00109 interpolate
00110 (
00111 const GeometricField<Type, fvPatchField, volMesh>& vf,
00112 const surfaceScalarField& faceFlux,
00113 Istream& schemeData
00114 )
00115 {
00116 if (surfaceInterpolation::debug)
00117 {
00118 Info<< "interpolate"
00119 << "(const GeometricField<Type, fvPatchField, volMesh>&, "
00120 << "const surfaceScalarField&, Istream&) : "
00121 << "interpolating GeometricField<Type, fvPatchField, volMesh> "
00122 << endl;
00123 }
00124
00125 return scheme<Type>(faceFlux, schemeData)().interpolate(vf);
00126 }
00127
00128
00129
00130 template<class Type>
00131 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
00132 interpolate
00133 (
00134 const GeometricField<Type, fvPatchField, volMesh>& vf,
00135 const surfaceScalarField& faceFlux,
00136 const word& name
00137 )
00138 {
00139 if (surfaceInterpolation::debug)
00140 {
00141 Info<< "interpolate"
00142 << "(const GeometricField<Type, fvPatchField, volMesh>&, "
00143 << "const surfaceScalarField&, const word&) : "
00144 << "interpolating GeometricField<Type, fvPatchField, volMesh> "
00145 << "using " << name
00146 << endl;
00147 }
00148
00149 return scheme<Type>(faceFlux, name)().interpolate(vf);
00150 }
00151
00152
00153 template<class Type>
00154 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
00155 interpolate
00156 (
00157 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf,
00158 const surfaceScalarField& faceFlux,
00159 const word& name
00160 )
00161 {
00162 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tsf =
00163 interpolate(tvf(), faceFlux, name);
00164
00165 tvf.clear();
00166
00167 return tsf;
00168 }
00169
00170
00171 template<class Type>
00172 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
00173 interpolate
00174 (
00175 const GeometricField<Type, fvPatchField, volMesh>& vf,
00176 const tmp<surfaceScalarField>& tFaceFlux,
00177 const word& name
00178 )
00179 {
00180 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tsf =
00181 interpolate(vf, tFaceFlux(), name);
00182
00183 tFaceFlux.clear();
00184
00185 return tsf;
00186 }
00187
00188
00189 template<class Type>
00190 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
00191 interpolate
00192 (
00193 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf,
00194 const tmp<surfaceScalarField>& tFaceFlux,
00195 const word& name
00196 )
00197 {
00198 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tsf =
00199 interpolate(tvf(), tFaceFlux(), name);
00200
00201 tvf.clear();
00202 tFaceFlux.clear();
00203
00204 return tsf;
00205 }
00206
00207
00208
00209 template<class Type>
00210 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
00211 interpolate
00212 (
00213 const GeometricField<Type, fvPatchField, volMesh>& vf,
00214 Istream& schemeData
00215 )
00216 {
00217 if (surfaceInterpolation::debug)
00218 {
00219 Info<< "interpolate"
00220 << "(const GeometricField<Type, fvPatchField, volMesh>&, "
00221 << "Istream&) : "
00222 << "interpolating GeometricField<Type, fvPatchField, volMesh> "
00223 << endl;
00224 }
00225
00226 return scheme<Type>(vf.mesh(), schemeData)().interpolate(vf);
00227 }
00228
00229
00230 template<class Type>
00231 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
00232 interpolate
00233 (
00234 const GeometricField<Type, fvPatchField, volMesh>& vf,
00235 const word& name
00236 )
00237 {
00238 if (surfaceInterpolation::debug)
00239 {
00240 Info<< "interpolate"
00241 << "(const GeometricField<Type, fvPatchField, volMesh>&, "
00242 << "const word&) : "
00243 << "interpolating GeometricField<Type, fvPatchField, volMesh> "
00244 << "using " << name
00245 << endl;
00246 }
00247
00248 return scheme<Type>(vf.mesh(), name)().interpolate(vf);
00249 }
00250
00251
00252 template<class Type>
00253 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
00254 interpolate
00255 (
00256 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf,
00257 const word& name
00258 )
00259 {
00260 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tsf =
00261 interpolate(tvf(), name);
00262
00263 tvf.clear();
00264
00265 return tsf;
00266 }
00267
00268
00269
00270 template<class Type>
00271 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
00272 interpolate
00273 (
00274 const GeometricField<Type, fvPatchField, volMesh>& vf
00275 )
00276 {
00277 if (surfaceInterpolation::debug)
00278 {
00279 Info<< "interpolate"
00280 << "(const GeometricField<Type, fvPatchField, volMesh>&) : "
00281 << "interpolating GeometricField<Type, fvPatchField, volMesh> "
00282 << "using run-time selected scheme"
00283 << endl;
00284 }
00285
00286 return interpolate(vf, "interpolate(" + vf.name() + ')');
00287 }
00288
00289
00290
00291 template<class Type>
00292 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
00293 interpolate
00294 (
00295 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf
00296 )
00297 {
00298 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tsf =
00299 interpolate(tvf());
00300 tvf.clear();
00301 return tsf;
00302 }
00303
00304
00305
00306
00307 }
00308
00309
00310
00311 }
00312
00313