Go to the documentation of this file.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
00027
00028
00029
00030
00031
00032
00033
00034
00035 #ifndef localBlended_H
00036 #define localBlended_H
00037
00038 #include <finiteVolume/surfaceInterpolationScheme.H>
00039
00040
00041
00042 namespace Foam
00043 {
00044
00045
00046
00047
00048
00049 template<class Type>
00050 class localBlended
00051 :
00052 public surfaceInterpolationScheme<Type>
00053 {
00054
00055
00056
00057 tmp<surfaceInterpolationScheme<Type> > tScheme1_;
00058
00059
00060 tmp<surfaceInterpolationScheme<Type> > tScheme2_;
00061
00062
00063
00064 localBlended(const localBlended&);
00065
00066
00067 void operator=(const localBlended&);
00068
00069
00070 public:
00071
00072
00073 TypeName("localBlended");
00074
00075
00076
00077
00078
00079
00080
00081 localBlended
00082 (
00083 const fvMesh& mesh,
00084 Istream& is
00085 )
00086 :
00087 surfaceInterpolationScheme<Type>(mesh),
00088 tScheme1_
00089 (
00090 surfaceInterpolationScheme<Type>::New(mesh, is)
00091 ),
00092 tScheme2_
00093 (
00094 surfaceInterpolationScheme<Type>::New(mesh, is)
00095 )
00096 {}
00097
00098
00099 localBlended
00100 (
00101 const fvMesh& mesh,
00102 const surfaceScalarField& faceFlux,
00103 Istream& is
00104 )
00105 :
00106 surfaceInterpolationScheme<Type>(mesh),
00107 tScheme1_
00108 (
00109 surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
00110 ),
00111 tScheme2_
00112 (
00113 surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
00114 )
00115 {}
00116
00117
00118
00119
00120
00121 tmp<surfaceScalarField> weights
00122 (
00123 const GeometricField<Type, fvPatchField, volMesh>& vf
00124 ) const
00125 {
00126 const surfaceScalarField& blendingFactor =
00127 this->mesh().objectRegistry::
lookupObject<const surfaceScalarField>
00128 (
00129 word(vf.name() + "BlendingFactor")
00130 );
00131
00132 return
00133 blendingFactor*tScheme1_().weights(vf)
00134 + (scalar(1) - blendingFactor)*tScheme2_().weights(vf);
00135 }
00136
00137
00138
00139 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
00140 interpolate(const GeometricField<Type, fvPatchField, volMesh>& vf) const
00141 {
00142 const surfaceScalarField& blendingFactor =
00143 (
00144 this->mesh().objectRegistry::
lookupObject<const surfaceScalarField>
00145 (
00146 word(vf.name() + "BlendingFactor")
00147 )
00148 );
00149
00150 return
00151 blendingFactor*tScheme1_().interpolate(vf)
00152 + (scalar(1) - blendingFactor)*tScheme2_().interpolate(vf);
00153 }
00154
00155
00156
00157 virtual bool corrected() const
00158 {
00159 return tScheme1_().corrected() || tScheme2_().corrected();
00160 }
00161
00162
00163
00164
00165 virtual tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
00166 correction
00167 (
00168 const GeometricField<Type, fvPatchField, volMesh>& vf
00169 ) const
00170 {
00171 const surfaceScalarField& blendingFactor =
00172 this->mesh().objectRegistry::
lookupObject<const surfaceScalarField>
00173 (
00174 word(vf.name() + "BlendingFactor")
00175 );
00176
00177 if (tScheme1_().corrected())
00178 {
00179 if (tScheme2_().corrected())
00180 {
00181 return
00182 (
00183 blendingFactor
00184 * tScheme1_().correction(vf)
00185 + (scalar(1.0) - blendingFactor)
00186 * tScheme2_().correction(vf)
00187 );
00188 }
00189 else
00190 {
00191 return
00192 (
00193 blendingFactor
00194 * tScheme1_().correction(vf)
00195 );
00196 }
00197 }
00198 else if (tScheme2_().corrected())
00199 {
00200 return
00201 (
00202 (scalar(1.0) - blendingFactor)
00203 * tScheme2_().correction(vf)
00204 );
00205 }
00206 else
00207 {
00208 return tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
00209 (
00210 NULL
00211 );
00212 }
00213 }
00214 };
00215
00216
00217
00218
00219 }
00220
00221
00222
00223 #endif
00224
00225
00226