mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: Added DimensionedField handling to volume/domain integrate
This commit is contained in:
@ -2,7 +2,7 @@
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
|
||||
\\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
@ -49,6 +49,7 @@ volumeIntegrate
|
||||
return vf.mesh().V()*vf.internalField();
|
||||
}
|
||||
|
||||
|
||||
template<class Type>
|
||||
tmp<Field<Type> >
|
||||
volumeIntegrate
|
||||
@ -62,6 +63,23 @@ volumeIntegrate
|
||||
}
|
||||
|
||||
|
||||
template<class Type>
|
||||
tmp<Field<Type> > volumeIntegrate(const DimensionedField<Type, volMesh>& df)
|
||||
{
|
||||
return df.mesh().V()*df.field();
|
||||
}
|
||||
|
||||
|
||||
template<class Type>
|
||||
tmp<Field<Type> >
|
||||
volumeIntegrate(const tmp<DimensionedField<Type, volMesh> >& tdf)
|
||||
{
|
||||
tmp<Field<Type> > tdidf = tdf().mesh().V()*tdf().field();
|
||||
tdf.clear();
|
||||
return tdidf;
|
||||
}
|
||||
|
||||
|
||||
template<class Type>
|
||||
dimensioned<Type>
|
||||
domainIntegrate
|
||||
@ -77,9 +95,9 @@ domainIntegrate
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
template<class Type>
|
||||
dimensioned<Type>
|
||||
domainIntegrate
|
||||
dimensioned<Type> domainIntegrate
|
||||
(
|
||||
const tmp<GeometricField<Type, fvPatchField, volMesh> >& tvf
|
||||
)
|
||||
@ -90,6 +108,33 @@ domainIntegrate
|
||||
}
|
||||
|
||||
|
||||
template<class Type>
|
||||
dimensioned<Type> domainIntegrate
|
||||
(
|
||||
const DimensionedField<Type, volMesh>& df
|
||||
)
|
||||
{
|
||||
return dimensioned<Type>
|
||||
(
|
||||
"domainIntegrate(" + df.name() + ')',
|
||||
dimVol*df.dimensions(),
|
||||
gSum(fvc::volumeIntegrate(df))
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
template<class Type>
|
||||
dimensioned<Type> domainIntegrate
|
||||
(
|
||||
const tmp<DimensionedField<Type, volMesh> >& tdf
|
||||
)
|
||||
{
|
||||
dimensioned<Type> integral = domainIntegrate(tdf());
|
||||
tdf.clear();
|
||||
return integral;
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace fvc
|
||||
|
||||
@ -2,7 +2,7 @@
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
|
||||
\\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
@ -67,6 +67,19 @@ namespace fvc
|
||||
);
|
||||
|
||||
|
||||
template<class Type>
|
||||
tmp<Field<Type> > volumeIntegrate
|
||||
(
|
||||
const DimensionedField<Type, volMesh>&
|
||||
);
|
||||
|
||||
template<class Type>
|
||||
tmp<Field<Type> > volumeIntegrate
|
||||
(
|
||||
const tmp<DimensionedField<Type, volMesh> >&
|
||||
);
|
||||
|
||||
|
||||
template<class Type>
|
||||
dimensioned<Type> domainIntegrate
|
||||
(
|
||||
@ -78,6 +91,19 @@ namespace fvc
|
||||
(
|
||||
const tmp<GeometricField<Type, fvPatchField, volMesh> >&
|
||||
);
|
||||
|
||||
|
||||
template<class Type>
|
||||
dimensioned<Type> domainIntegrate
|
||||
(
|
||||
const DimensionedField<Type, volMesh>&
|
||||
);
|
||||
|
||||
template<class Type>
|
||||
dimensioned<Type> domainIntegrate
|
||||
(
|
||||
const tmp<DimensionedField<Type, volMesh> >&
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
|
||||
Reference in New Issue
Block a user