Updates to the fieldValues function object

- Updates to enable correct operation in parallel
- Added weighted average operation for cell sources
This commit is contained in:
andy
2010-01-08 17:14:53 +00:00
parent 711e0a42dd
commit 1812773b25
8 changed files with 249 additions and 174 deletions

View File

@ -220,7 +220,7 @@ void Foam::fieldValues::faceSource::initialise(const dictionary& dict)
}
default:
{
FatalErrorIn("faceSource::constructFaceAddressing()")
FatalErrorIn("faceSource::initiliase()")
<< type() << " " << name_ << ": "
<< sourceTypeNames_[source_] << "(" << sourceName_ << "):"
<< nl << " Unknown source type. Valid source types are:"
@ -245,7 +245,7 @@ void Foam::fieldValues::faceSource::initialise(const dictionary& dict)
}
else
{
FatalErrorIn("faceSource::constructFaceAddressing()")
FatalErrorIn("faceSource::initialise()")
<< type() << " " << name_ << ": "
<< sourceTypeNames_[source_] << "(" << sourceName_ << "):"
<< nl << " Weight field " << weightFieldName_
@ -326,9 +326,12 @@ void Foam::fieldValues::faceSource::write()
if (active_)
{
outputFilePtr_()
<< obr_.time().value() << tab
<< sum(filterField(mesh().magSf()));
if (Pstream::master())
{
outputFilePtr_()
<< obr_.time().value() << tab
<< sum(filterField(mesh().magSf()));
}
forAll(fields_, i)
{
@ -339,7 +342,10 @@ void Foam::fieldValues::faceSource::write()
writeValues<tensor>(fields_[i]);
}
outputFilePtr_()<< endl;
if (Pstream::master())
{
outputFilePtr_()<< endl;
}
if (log_)
{
@ -350,4 +356,3 @@ void Foam::fieldValues::faceSource::write()
// ************************************************************************* //

View File

@ -153,17 +153,22 @@ protected:
//- Initialise, e.g. face addressing
void initialise(const dictionary& dict);
//- Insert field values into values list
//- Return true if the field name is valid
template<class Type>
bool setFieldValues
(
const word& fieldName,
List<Type>& values
) const;
bool validField(const word& fieldName) const;
//- Return field values by looking up field name
template<class Type>
tmp<Field<Type> > setFieldValues(const word& fieldName) const;
//- Apply the 'operation' to the values
template<class Type>
Type processValues(const List<Type>& values) const;
Type processValues
(
const Field<Type>& values,
const scalarField& magSf,
const scalarField& weightField
) const;
//- Output file header information
virtual void writeFileHeader();

View File

@ -33,70 +33,17 @@ License
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
template<class Type>
bool Foam::fieldValues::faceSource::setFieldValues
(
const word& fieldName,
List<Type>& values
) const
bool Foam::fieldValues::faceSource::validField(const word& fieldName) const
{
values.setSize(faceId_.size(), pTraits<Type>::zero);
typedef GeometricField<Type, fvsPatchField, surfaceMesh> sf;
typedef GeometricField<Type, fvPatchField, volMesh> vf;
if (obr_.foundObject<sf>(fieldName))
{
const sf& field = obr_.lookupObject<sf>(fieldName);
forAll(values, i)
{
label faceI = faceId_[i];
label patchI = facePatchId_[i];
if (patchI >= 0)
{
values[i] = field.boundaryField()[patchI][faceI];
}
else
{
values[i] = field[faceI];
}
values[i] *= flipMap_[i];
}
return true;
}
else if (obr_.foundObject<vf>(fieldName))
{
const vf& field = obr_.lookupObject<vf>(fieldName);
forAll(values, i)
{
label faceI = faceId_[i];
label patchI = facePatchId_[i];
if (patchI >= 0)
{
values[i] = field.boundaryField()[patchI][faceI];
}
else
{
FatalErrorIn
(
"fieldValues::faceSource::setFieldValues"
"("
"const word&, "
"List<Type>&"
") const"
) << type() << " " << name_ << ": "
<< sourceTypeNames_[source_] << "(" << sourceName_ << "):"
<< nl
<< " Unable to process internal faces for volume field "
<< fieldName << nl << abort(FatalError);
}
values[i] *= flipMap_[i];
}
return true;
}
@ -104,10 +51,34 @@ bool Foam::fieldValues::faceSource::setFieldValues
}
template<class Type>
Foam::tmp<Foam::Field<Type> > Foam::fieldValues::faceSource::setFieldValues
(
const word& fieldName
) const
{
typedef GeometricField<Type, fvsPatchField, surfaceMesh> sf;
typedef GeometricField<Type, fvPatchField, volMesh> vf;
if (obr_.foundObject<sf>(fieldName))
{
return filterField(obr_.lookupObject<sf>(fieldName));
}
else if (obr_.foundObject<vf>(fieldName))
{
return filterField(obr_.lookupObject<vf>(fieldName));
}
return tmp<Field<Type> >(new Field<Type>(0.0));
}
template<class Type>
Type Foam::fieldValues::faceSource::processValues
(
const List<Type>& values
const Field<Type>& values,
const scalarField& magSf,
const scalarField& weightField
) const
{
Type result = pTraits<Type>::zero;
@ -120,54 +91,17 @@ Type Foam::fieldValues::faceSource::processValues
}
case opAreaAverage:
{
tmp<scalarField> magSf = filterField(mesh().magSf());
result = sum(values*magSf())/sum(magSf());
result = sum(values*magSf)/sum(magSf);
break;
}
case opAreaIntegrate:
{
result = sum(values*filterField(mesh().magSf()));
result = sum(values*magSf);
break;
}
case opWeightedAverage:
{
if (mesh().foundObject<volScalarField>(weightFieldName_))
{
tmp<scalarField> wField =
filterField
(
mesh().lookupObject<volScalarField>(weightFieldName_)
);
result = sum(values*wField())/sum(wField());
}
else if (mesh().foundObject<surfaceScalarField>(weightFieldName_))
{
tmp<scalarField> wField =
filterField
(
mesh().lookupObject<surfaceScalarField>
(
weightFieldName_
)
);
result = sum(values*wField())/sum(wField());
}
else
{
FatalErrorIn
(
"fieldValues::faceSource::processValues"
"("
"List<Type>&"
") const"
) << type() << " " << name_ << ": "
<< sourceTypeNames_[source_] << "(" << sourceName_ << "):"
<< nl
<< " Weight field " << weightFieldName_
<< " must be either a " << volScalarField::typeName
<< " or " << surfaceScalarField::typeName << nl
<< abort(FatalError);
}
result = sum(values*weightField)/sum(weightField);
break;
}
default:
@ -185,25 +119,20 @@ Type Foam::fieldValues::faceSource::processValues
template<class Type>
bool Foam::fieldValues::faceSource::writeValues(const word& fieldName)
{
List<List<Type> > allValues(Pstream::nProcs());
const bool ok = validField<Type>(fieldName);
bool validField =
setFieldValues<Type>(fieldName, allValues[Pstream::myProcNo()]);
if (validField)
if (ok)
{
Pstream::gatherList(allValues);
Field<Type> values = combineFields(setFieldValues<Type>(fieldName));
scalarField magSf = combineFields(filterField(mesh().magSf()));
scalarField weightField =
combineFields(setFieldValues<scalar>(weightFieldName_));
if (Pstream::master())
{
List<Type> values =
ListListOps::combine<List<Type> >
(
allValues,
accessOp<List<Type> >()
);
Type result = processValues(values);
Type result = processValues(values, magSf, weightField);
if (valueOutput_)
{
@ -222,7 +151,6 @@ bool Foam::fieldValues::faceSource::writeValues(const word& fieldName)
).write();
}
outputFilePtr_()<< tab << result;
if (log_)
@ -234,7 +162,7 @@ bool Foam::fieldValues::faceSource::writeValues(const word& fieldName)
}
}
return validField;
return ok;
}