/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2011-2013 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 "processorCyclicPointPatchField.H" #include "transformField.H" #include "processorPolyPatch.H" // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // template Foam::processorCyclicPointPatchField::processorCyclicPointPatchField ( const pointPatch& p, const DimensionedField& iF ) : coupledPointPatchField(p, iF), procPatch_(refCast(p)), receiveBuf_(0) {} template Foam::processorCyclicPointPatchField::processorCyclicPointPatchField ( const pointPatch& p, const DimensionedField& iF, const dictionary& dict ) : coupledPointPatchField(p, iF, dict), procPatch_(refCast(p)), receiveBuf_(0) {} template Foam::processorCyclicPointPatchField::processorCyclicPointPatchField ( const processorCyclicPointPatchField& ptf, const pointPatch& p, const DimensionedField& iF, const pointPatchFieldMapper& mapper ) : coupledPointPatchField(ptf, p, iF, mapper), procPatch_(refCast(ptf.patch())), receiveBuf_(0) {} template Foam::processorCyclicPointPatchField::processorCyclicPointPatchField ( const processorCyclicPointPatchField& ptf, const DimensionedField& iF ) : coupledPointPatchField(ptf, iF), procPatch_(refCast(ptf.patch())), receiveBuf_(0) {} // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * // template Foam::processorCyclicPointPatchField::~processorCyclicPointPatchField() {} // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // template void Foam::processorCyclicPointPatchField::initSwapAddSeparated ( const Pstream::commsTypes commsType, Field& pField ) const { if (Pstream::parRun()) { // Get internal field into correct order for opposite side Field pf ( this->patchInternalField ( pField, procPatch_.reverseMeshPoints() ) ); if (commsType == Pstream::nonBlocking) { receiveBuf_.setSize(pf.size()); IPstream::read ( commsType, procPatch_.neighbProcNo(), reinterpret_cast(receiveBuf_.begin()), receiveBuf_.byteSize(), procPatch_.tag(), procPatch_.comm() ); } OPstream::write ( commsType, procPatch_.neighbProcNo(), reinterpret_cast(pf.begin()), pf.byteSize(), procPatch_.tag(), procPatch_.comm() ); } } template void Foam::processorCyclicPointPatchField::swapAddSeparated ( const Pstream::commsTypes commsType, Field& pField ) const { if (Pstream::parRun()) { // If nonblocking data has already been received into receiveBuf_ if (commsType != Pstream::nonBlocking) { receiveBuf_.setSize(this->size()); IPstream::read ( commsType, procPatch_.neighbProcNo(), reinterpret_cast(receiveBuf_.begin()), receiveBuf_.byteSize(), procPatch_.tag(), procPatch_.comm() ); } if (doTransform()) { const processorCyclicPolyPatch& ppp = procPatch_.procCyclicPolyPatch(); const tensor& forwardT = ppp.forwardT()[0]; transform(receiveBuf_, forwardT, receiveBuf_); } // All points are separated this->addToInternalField(pField, receiveBuf_); } } // ************************************************************************* //