/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2012-2016 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License This file is part of OpenFOAM. OpenFOAM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OpenFOAM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OpenFOAM. If not, see . \*---------------------------------------------------------------------------*/ #include "fvMesh.H" #include "volFields.H" #include "directFvPatchFieldMapper.H" #include "calculatedFvPatchField.H" #include "weightedFvPatchFieldMapper.H" // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // namespace Foam { //- Helper class for list template class ListPlusEqOp { public: void operator()(List& x, const List y) const { if (y.size()) { if (x.size()) { label sz = x.size(); x.setSize(sz + y.size()); forAll(y, i) { x[sz++] = y[i]; } } else { x = y; } } } }; } // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * // template void Foam::meshToMesh::add ( UList& fld, const label offset ) const { forAll(fld, i) { fld[i] += offset; } } template void Foam::meshToMesh::mapSrcToTgt ( const UList& srcField, const CombineOp& cop, List& result ) const { if (result.size() != tgtToSrcCellAddr_.size()) { FatalErrorInFunction << "Supplied field size is not equal to target mesh size" << nl << " source mesh = " << srcToTgtCellAddr_.size() << nl << " target mesh = " << tgtToSrcCellAddr_.size() << nl << " supplied field = " << result.size() << abort(FatalError); } multiplyWeightedOp cbop(cop); if (singleMeshProc_ == -1) { const mapDistribute& map = srcMapPtr_(); List work(srcField); map.distribute(work); forAll(result, cellI) { const labelList& srcAddress = tgtToSrcCellAddr_[cellI]; const scalarList& srcWeight = tgtToSrcCellWght_[cellI]; if (srcAddress.size()) { // result[cellI] = Zero; result[cellI] *= (1.0 - sum(srcWeight)); forAll(srcAddress, i) { label srcI = srcAddress[i]; scalar w = srcWeight[i]; cbop(result[cellI], cellI, work[srcI], w); } } } } else { forAll(result, cellI) { const labelList& srcAddress = tgtToSrcCellAddr_[cellI]; const scalarList& srcWeight = tgtToSrcCellWght_[cellI]; if (srcAddress.size()) { // result[cellI] = Zero; result[cellI] *= (1.0 - sum(srcWeight)); forAll(srcAddress, i) { label srcI = srcAddress[i]; scalar w = srcWeight[i]; cbop(result[cellI], cellI, srcField[srcI], w); } } } } } template Foam::tmp> Foam::meshToMesh::mapSrcToTgt ( const Field& srcField, const CombineOp& cop ) const { tmp> tresult ( new Field ( tgtToSrcCellAddr_.size(), Zero ) ); mapSrcToTgt(srcField, cop, tresult.ref()); return tresult; } template Foam::tmp> Foam::meshToMesh::mapSrcToTgt ( const tmp>& tsrcField, const CombineOp& cop ) const { return mapSrcToTgt(tsrcField(), cop); } template Foam::tmp> Foam::meshToMesh::mapSrcToTgt ( const Field& srcField ) const { return mapSrcToTgt(srcField, plusEqOp()); } template Foam::tmp> Foam::meshToMesh::mapSrcToTgt ( const tmp>& tsrcField ) const { return mapSrcToTgt(tsrcField()); } template void Foam::meshToMesh::mapTgtToSrc ( const UList& tgtField, const CombineOp& cop, List& result ) const { if (result.size() != srcToTgtCellAddr_.size()) { FatalErrorInFunction << "Supplied field size is not equal to source mesh size" << nl << " source mesh = " << srcToTgtCellAddr_.size() << nl << " target mesh = " << tgtToSrcCellAddr_.size() << nl << " supplied field = " << result.size() << abort(FatalError); } multiplyWeightedOp cbop(cop); if (singleMeshProc_ == -1) { const mapDistribute& map = tgtMapPtr_(); List work(tgtField); map.distribute(work); forAll(result, cellI) { const labelList& tgtAddress = srcToTgtCellAddr_[cellI]; const scalarList& tgtWeight = srcToTgtCellWght_[cellI]; if (tgtAddress.size()) { result[cellI] *= (1.0 - sum(tgtWeight)); forAll(tgtAddress, i) { label tgtI = tgtAddress[i]; scalar w = tgtWeight[i]; cbop(result[cellI], cellI, work[tgtI], w); } } } } else { forAll(result, cellI) { const labelList& tgtAddress = srcToTgtCellAddr_[cellI]; const scalarList& tgtWeight = srcToTgtCellWght_[cellI]; if (tgtAddress.size()) { result[cellI] *= (1.0 - sum(tgtWeight)); forAll(tgtAddress, i) { label tgtI = tgtAddress[i]; scalar w = tgtWeight[i]; cbop(result[cellI], cellI, tgtField[tgtI], w); } } } } } template Foam::tmp> Foam::meshToMesh::mapTgtToSrc ( const Field& tgtField, const CombineOp& cop ) const { tmp> tresult ( new Field ( srcToTgtCellAddr_.size(), Zero ) ); mapTgtToSrc(tgtField, cop, tresult.ref()); return tresult; } template Foam::tmp> Foam::meshToMesh::mapTgtToSrc ( const tmp>& ttgtField, const CombineOp& cop ) const { return mapTgtToSrc(ttgtField(), cop); } template Foam::tmp> Foam::meshToMesh::mapTgtToSrc ( const Field& tgtField ) const { return mapTgtToSrc(tgtField, plusEqOp()); } template Foam::tmp> Foam::meshToMesh::mapTgtToSrc ( const tmp>& ttgtField ) const { return mapTgtToSrc(ttgtField(), plusEqOp()); } template void Foam::meshToMesh::mapSrcToTgt ( const GeometricField& field, const CombineOp& cop, GeometricField& result ) const { mapSrcToTgt(field, cop, result.internalField()); const PtrList& AMIList = patchAMIs(); forAll(AMIList, i) { label srcPatchI = srcPatchID_[i]; label tgtPatchI = tgtPatchID_[i]; const fvPatchField& srcField = field.boundaryField()[srcPatchI]; fvPatchField& tgtField = result.boundaryField()[tgtPatchI]; // 2.3 does not do distributed mapping yet so only do if // running on single processor if (AMIList[i].singlePatchProc() != -1) { // Clone and map (since rmap does not do general mapping) tmp> tnewTgt ( fvPatchField::New ( srcField, tgtField.patch(), result.dimensionedInternalField(), weightedFvPatchFieldMapper ( AMIList[i].tgtAddress(), AMIList[i].tgtWeights() ) ) ); // Transfer all mapped quantities (value and e.g. gradient) onto // tgtField. Value will get overwritten below. tgtField.rmap(tnewTgt(), identity(tgtField.size())); } tgtField == Type(Zero); AMIList[i].interpolateToTarget ( srcField, multiplyWeightedOp(cop), tgtField, UList::null() ); } forAll(cuttingPatches_, i) { label patchI = cuttingPatches_[i]; fvPatchField& pf = result.boundaryField()[patchI]; pf == pf.patchInternalField(); } } template Foam::tmp> Foam::meshToMesh::mapSrcToTgt ( const GeometricField& field, const CombineOp& cop ) const { typedef GeometricField fieldType; const fvMesh& tgtMesh = static_cast(tgtRegion_); const fvBoundaryMesh& tgtBm = tgtMesh.boundary(); const typename fieldType::GeometricBoundaryField& srcBfld = field.boundaryField(); PtrList> tgtPatchFields(tgtBm.size()); // constuct tgt boundary patch types as copy of 'field' boundary types // note: this will provide place holders for fields with additional // entries, but these values will need to be reset forAll(tgtPatchID_, i) { label srcPatchI = srcPatchID_[i]; label tgtPatchI = tgtPatchID_[i]; if (!tgtPatchFields.set(tgtPatchI)) { tgtPatchFields.set ( tgtPatchI, fvPatchField::New ( srcBfld[srcPatchI], tgtMesh.boundary()[tgtPatchI], DimensionedField::null(), directFvPatchFieldMapper ( labelList(tgtMesh.boundary()[tgtPatchI].size(), -1) ) ) ); } } // Any unset tgtPatchFields become calculated forAll(tgtPatchFields, tgtPatchI) { if (!tgtPatchFields.set(tgtPatchI)) { // Note: use factory New method instead of direct generation of // calculated so we keep constraints tgtPatchFields.set ( tgtPatchI, fvPatchField::New ( calculatedFvPatchField::typeName, tgtMesh.boundary()[tgtPatchI], DimensionedField::null() ) ); } } tmp tresult ( new fieldType ( IOobject ( type() + ":interpolate(" + field.name() + ")", tgtMesh.time().timeName(), tgtMesh, IOobject::NO_READ, IOobject::NO_WRITE ), tgtMesh, field.dimensions(), Field(tgtMesh.nCells(), Zero), tgtPatchFields ) ); mapSrcToTgt(field, cop, tresult.ref()); return tresult; } template Foam::tmp> Foam::meshToMesh::mapSrcToTgt ( const tmp>& tfield, const CombineOp& cop ) const { return mapSrcToTgt(tfield(), cop); } template Foam::tmp> Foam::meshToMesh::mapSrcToTgt ( const GeometricField& field ) const { return mapSrcToTgt(field, plusEqOp()); } template Foam::tmp> Foam::meshToMesh::mapSrcToTgt ( const tmp>& tfield ) const { return mapSrcToTgt(tfield(), plusEqOp()); } template void Foam::meshToMesh::mapTgtToSrc ( const GeometricField& field, const CombineOp& cop, GeometricField& result ) const { mapTgtToSrc(field, cop, result.internalField()); const PtrList& AMIList = patchAMIs(); forAll(AMIList, i) { label srcPatchI = srcPatchID_[i]; label tgtPatchI = tgtPatchID_[i]; fvPatchField& srcField = result.boundaryField()[srcPatchI]; const fvPatchField& tgtField = field.boundaryField()[tgtPatchI]; // 2.3 does not do distributed mapping yet so only do if // running on single processor if (AMIList[i].singlePatchProc() != -1) { // Clone and map (since rmap does not do general mapping) tmp> tnewSrc ( fvPatchField::New ( tgtField, srcField.patch(), result.dimensionedInternalField(), weightedFvPatchFieldMapper ( AMIList[i].srcAddress(), AMIList[i].srcWeights() ) ) ); // Transfer all mapped quantities (value and e.g. gradient) onto // srcField. Value will get overwritten below srcField.rmap(tnewSrc(), identity(srcField.size())); } srcField == Type(Zero); AMIList[i].interpolateToSource ( tgtField, multiplyWeightedOp(cop), srcField, UList::null() ); } forAll(cuttingPatches_, i) { label patchI = cuttingPatches_[i]; fvPatchField& pf = result.boundaryField()[patchI]; pf == pf.patchInternalField(); } } template Foam::tmp> Foam::meshToMesh::mapTgtToSrc ( const GeometricField& field, const CombineOp& cop ) const { typedef GeometricField fieldType; const fvMesh& srcMesh = static_cast(srcRegion_); const fvBoundaryMesh& srcBm = srcMesh.boundary(); const typename fieldType::GeometricBoundaryField& tgtBfld = field.boundaryField(); PtrList> srcPatchFields(srcBm.size()); // constuct src boundary patch types as copy of 'field' boundary types // note: this will provide place holders for fields with additional // entries, but these values will need to be reset forAll(srcPatchID_, i) { label srcPatchI = srcPatchID_[i]; label tgtPatchI = tgtPatchID_[i]; if (!srcPatchFields.set(tgtPatchI)) { srcPatchFields.set ( srcPatchI, fvPatchField::New ( tgtBfld[srcPatchI], srcMesh.boundary()[tgtPatchI], DimensionedField::null(), directFvPatchFieldMapper ( labelList(srcMesh.boundary()[srcPatchI].size(), -1) ) ) ); } } // Any unset srcPatchFields become calculated forAll(srcPatchFields, srcPatchI) { if (!srcPatchFields.set(srcPatchI)) { // Note: use factory New method instead of direct generation of // calculated so we keep constraints srcPatchFields.set ( srcPatchI, fvPatchField::New ( calculatedFvPatchField::typeName, srcMesh.boundary()[srcPatchI], DimensionedField::null() ) ); } } tmp tresult ( new fieldType ( IOobject ( type() + ":interpolate(" + field.name() + ")", srcMesh.time().timeName(), srcMesh, IOobject::NO_READ, IOobject::NO_WRITE ), srcMesh, field.dimensions(), Field(srcMesh.nCells(), Zero), srcPatchFields ) ); mapTgtToSrc(field, cop, tresult.ref()); return tresult; } template Foam::tmp> Foam::meshToMesh::mapTgtToSrc ( const tmp>& tfield, const CombineOp& cop ) const { return mapTgtToSrc(tfield(), cop); } template Foam::tmp> Foam::meshToMesh::mapTgtToSrc ( const GeometricField& field ) const { return mapTgtToSrc(field, plusEqOp()); } template Foam::tmp> Foam::meshToMesh::mapTgtToSrc ( const tmp>& tfield ) const { return mapTgtToSrc(tfield(), plusEqOp()); } // ************************************************************************* //