Report patch name instead of index in debug
[foam-extend-3.2.git] / src / foam / meshes / polyMesh / polyPatches / constraint / ggi / ggiPolyPatch.C
blob99ee1b8c53ad3b13cd17ddfe9658ced959583278
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation, either version 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Author
25     Hrvoje Jasak, Wikki Ltd.  All rights reserved.
27 Contributor
28     Martin Beaudoin, Hydro-Quebec, (2008)
30 \*---------------------------------------------------------------------------*/
32 #include "ggiPolyPatch.H"
33 #include "polyBoundaryMesh.H"
34 #include "addToRunTimeSelectionTable.H"
35 #include "demandDrivenData.H"
36 #include "polyPatchID.H"
37 #include "ZoneIDs.H"
38 #include "SubField.H"
39 #include "foamTime.H"
40 #include "indirectPrimitivePatch.H"
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 namespace Foam
46     defineTypeNameAndDebugWithDescription
47     (
48         ggiPolyPatch,
49         0,
50         "If value > 1, write uncovered GGI patch facets to VTK file"
51     );
53     addToRunTimeSelectionTable(polyPatch, ggiPolyPatch, word);
54     addToRunTimeSelectionTable(polyPatch, ggiPolyPatch, dictionary);
58 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
60 bool Foam::ggiPolyPatch::active() const
62     polyPatchID shadow(shadowName_, boundaryMesh());
63     faceZoneID zone(zoneName_, boundaryMesh().mesh().faceZones());
65     return shadow.active() && zone.active();
69 void Foam::ggiPolyPatch::calcZoneAddressing() const
71     // Calculate patch-to-zone addressing
72     if (zoneAddressingPtr_)
73     {
74         FatalErrorIn("void ggiPolyPatch::calcZoneAddressing() const")
75             << "Patch to zone addressing already calculated"
76             << abort(FatalError);
77     }
79     // Calculate patch-to-zone addressing
80     zoneAddressingPtr_ = new labelList(size());
81     labelList& zAddr = *zoneAddressingPtr_;
82     const faceZone& myZone = zone();
84     for (label i = 0; i < size(); i++)
85     {
86         zAddr[i] = myZone.whichFace(start() + i);
87     }
89     // Check zone addressing
90     if (zAddr.size() > 0 && min(zAddr) < 0)
91     {
92         Info<< "myZone: " << myZone << nl
93             << "my start and size: " << start() << " and " << size() << nl
94             << "zAddr: " << zAddr << endl;
96         FatalErrorIn("void ggiPolyPatch::calcZoneAddressing() const")
97             << "Problem with patch-to-zone addressing: some patch faces "
98             << "not found in interpolation zone"
99             << abort(FatalError);
100     }
104 void Foam::ggiPolyPatch::calcRemoteZoneAddressing() const
106     // Calculate patch-to-zone addressing
107     if (remoteZoneAddressingPtr_)
108     {
109         FatalErrorIn("void ggiPolyPatch::calcRemoteZoneAddressing() const")
110             << "Patch to remote zone addressing already calculated"
111             << abort(FatalError);
112     }
114     if (debug)
115     {
116         Pout<< "ggiPolyPatch::calcRemoteZoneAddressing() const for patch "
117             << name() << endl;
118     }
120     // Once zone addressing is established, visit the opposite side and find
121     // out which face data is needed for interpolation
122     boolList usedShadows(shadow().zone().size(), false);
124     const labelList& zAddr = zoneAddressing();
126     if (master())
127     {
128         const labelListList& addr = patchToPatch().masterAddr();
130         forAll (zAddr, mfI)
131         {
132             const labelList& nbrs = addr[zAddr[mfI]];
134             forAll (nbrs, nbrI)
135             {
136                 usedShadows[nbrs[nbrI]] = true;
137             }
138         }
139     }
140     else
141     {
142         const labelListList& addr = patchToPatch().slaveAddr();
144         forAll (zAddr, mfI)
145         {
146             const labelList& nbrs = addr[zAddr[mfI]];
148             forAll (nbrs, nbrI)
149             {
150                 usedShadows[nbrs[nbrI]] = true;
151             }
152         }
153     }
155     // Count and pick up shadow indices
156     label nShadows = 0;
158     forAll (usedShadows, sI)
159     {
160         if (usedShadows[sI])
161         {
162             nShadows++;
163         }
164     }
166     remoteZoneAddressingPtr_ = new labelList(nShadows);
167     labelList& rza = *remoteZoneAddressingPtr_;
169     // Reset counter for re-use
170     nShadows = 0;
172     forAll (usedShadows, sI)
173     {
174         if (usedShadows[sI])
175         {
176             rza[nShadows] = sI;
177             nShadows++;
178         }
179     }
183 void Foam::ggiPolyPatch::calcPatchToPatch() const
185     // Create patch-to-patch interpolation
186     if (patchToPatchPtr_)
187     {
188         FatalErrorIn("void ggiPolyPatch::calcPatchToPatch() const")
189             << "Patch to patch interpolation already calculated"
190             << abort(FatalError);
191     }
193     if (master())
194     {
195         if (debug)
196         {
197             InfoIn("void ggiPolyPatch::calcPatchToPatch() const")
198                 << "Calculating patch to patch interpolation for patch"
199                 << name() << endl;
200         }
202         // Create interpolation for zones
203         patchToPatchPtr_ =
204             new ggiZoneInterpolation
205             (
206                 zone()(),           // This zone reference
207                 shadow().zone()(),  // Shadow zone reference
208                 forwardT(),
209                 reverseT(),
210                 -separation(), // Slave-to-master separation: Use - localValue
211                 // Bug fix, delayed slave evaluation causes error
212                 // HJ, 30/Jun/2013
213                 0,             // Non-overlapping face tolerances
214                 0,             // HJ, 24/Oct/2008
215                 true,          // Rescale weighting factors.  Bug fix, MB.
216                 reject_        // Quick rejection algorithm, default BB_OCTREE
217             );
219         // Abort immediately if uncovered faces are present and the option
220         // bridgeOverlap is not set.
221         if
222         (
223             (
224                 patchToPatch().uncoveredMasterFaces().size() > 0
225             && !bridgeOverlap()
226             )
227          || (
228                 patchToPatch().uncoveredSlaveFaces().size() > 0
229             && !shadow().bridgeOverlap()
230             )
231         )
232         {
233             FatalErrorIn("void ggiPolyPatch::calcPatchToPatch() const")
234                 << "Found uncovered faces for GGI interface "
235                 << name() << "/" << shadowName()
236                 << " while the bridgeOverlap option is not set "
237                 << "in the boundary file." << endl
238                 << "This is an unrecoverable error. Aborting."
239                 << abort(FatalError);
240         }
241     }
242     else
243     {
244         FatalErrorIn("void ggiPolyPatch::calcPatchToPatch() const")
245             << "Attempting to create GGIInterpolation on a shadow"
246             << abort(FatalError);
247     }
251 void Foam::ggiPolyPatch::calcReconFaceCellCentres() const
253     if (reconFaceCellCentresPtr_)
254     {
255         FatalErrorIn
256         (
257             "void ggiPolyPatch::calcReconFaceCellCentres() const"
258         )   << "Reconstructed cell centres already calculated"
259             << abort(FatalError);
260     }
262     if (debug)
263     {
264         InfoIn("void ggiPolyPatch::calcReconFaceCellCentres() const")
265             << "Calculating recon centres for patch"
266             << name() << endl;
267     }
269     // Create neighbouring face centres using interpolation
270     if (master())
271     {
272         const label shadowID = shadowIndex();
274         // Get the transformed and interpolated shadow face cell centers
275         reconFaceCellCentresPtr_ =
276             new vectorField
277             (
278                 interpolate
279                 (
280                     boundaryMesh()[shadowID].faceCellCentres()
281                   - boundaryMesh()[shadowID].faceCentres()
282                 )
283               + faceCentres()
284             );
285     }
286     else
287     {
288         FatalErrorIn("void ggiPolyPatch::calcReconFaceCellCentres() const")
289             << "Attempting to create reconFaceCellCentres on a shadow"
290             << abort(FatalError);
291     }
295 void Foam::ggiPolyPatch::calcLocalParallel() const
297     // Calculate patch-to-zone addressing
298     if (localParallelPtr_)
299     {
300         FatalErrorIn("void ggiPolyPatch::calcLocalParallel() const")
301             << "Local parallel switch already calculated"
302             << abort(FatalError);
303     }
305     localParallelPtr_ = new bool(false);
306     bool& emptyOrComplete = *localParallelPtr_;
308     // If running in serial, all GGIs are expanded to zone size.
309     // This happens on decomposition and reconstruction where
310     // size and shadow size may be zero, but zone size may not
311     // HJ, 1/Jun/2011
312     if (!Pstream::parRun())
313     {
314         emptyOrComplete = false;
315     }
316     else
317     {
318         // Check that patch size is greater than the zone size.
319         // This is an indication of the error where the face zone is not global
320         // in a parallel run.  HJ, 9/Nov/2014
321         if (size() > zone().size())
322         {
323             FatalErrorIn("void ggiPolyPatch::calcLocalParallel() const")
324                 << "Patch size is greater than zone size for GGI patch "
325                 << name() << ".  This is not allowerd: "
326                 << "the face zone must contain all patch faces and be "
327                 << "global in parallel runs"
328                 << abort(FatalError);
329         }
331         // Calculate localisation on master and shadow
332         emptyOrComplete =
333             (
334                 zone().size() == size()
335              && shadow().zone().size() == shadow().size()
336             )
337          || (size() == 0 && shadow().size() == 0);
339         reduce(emptyOrComplete, andOp<bool>());
340     }
342     if (debug && Pstream::parRun())
343     {
344         Info<< "GGI patch Master: " << name()
345             << " Slave: " << shadowName() << " is ";
347         if (emptyOrComplete)
348         {
349            Info<< "local parallel" << endl;
350         }
351         else
352         {
353             Info<< "split between multiple processors" << endl;
354         }
355     }
359 void Foam::ggiPolyPatch::calcSendReceive() const
361     // Note: all processors will execute calcSendReceive but only master will
362     // hold the information.  Therefore, pointers on slave processors
363     // will remain meaningless, but for purposes of consistency
364     // (of the calc-call) they will be set to zero-sized array
365     // HJ, 4/Jun/2011
367     if (receiveAddrPtr_ || sendAddrPtr_)
368     {
369         FatalErrorIn("void ggiPolyPatch::calcSendReceive() const")
370             << "Send-receive addressing already calculated"
371             << abort(FatalError);
372     }
374     if (debug)
375     {
376         Pout<< "ggiPolyPatch::calcSendReceive() const for patch "
377             << name() << endl;
378     }
380     if (!Pstream::parRun())
381     {
382         FatalErrorIn("void ggiPolyPatch::calcSendReceive() const")
383             << "Requested calculation of send-receive addressing for a "
384             << "serial run.  This is not allowed"
385             << abort(FatalError);
386     }
388     // Master will receive and store the maps
389     if (Pstream::master())
390     {
391         receiveAddrPtr_ = new labelListList(Pstream::nProcs());
392         labelListList& rAddr = *receiveAddrPtr_;
394         sendAddrPtr_ = new labelListList(Pstream::nProcs());
395         labelListList& sAddr = *sendAddrPtr_;
397         // Insert master
398         rAddr[0] = zoneAddressing();
400         for (label procI = 1; procI < Pstream::nProcs(); procI++)
401         {
402             // Note: must use normal comms because the size of the
403             // communicated lists is unknown on the receiving side
404             // HJ, 4/Jun/2011
406             // Opt: reconsider mode of communication
407             IPstream ip(Pstream::scheduled, procI);
409             rAddr[procI] = labelList(ip);
411             sAddr[procI] = labelList(ip);
412         }
413     }
414     else
415     {
416         // Create dummy pointers: only master processor stores maps
417         receiveAddrPtr_ = new labelListList();
418         sendAddrPtr_ = new labelListList();
420         // Send information to master
421         const labelList& za = zoneAddressing();
422         const labelList& ra = remoteZoneAddressing();
424         // Note: must use normal comms because the size of the
425         // communicated lists is unknown on the receiving side
426         // HJ, 4/Jun/2011
428         // Opt: reconsider mode of communication
429         OPstream op(Pstream::scheduled, Pstream::masterNo());
431         // Send local and remote addressing to master
432         op << za << ra;
433     }
437 const Foam::labelListList& Foam::ggiPolyPatch::receiveAddr() const
439     if (!receiveAddrPtr_)
440     {
441         calcSendReceive();
442     }
444     return *receiveAddrPtr_;
448 const Foam::labelListList& Foam::ggiPolyPatch::sendAddr() const
450     if (!sendAddrPtr_)
451     {
452         calcSendReceive();
453     }
455     return *sendAddrPtr_;
459 void Foam::ggiPolyPatch::clearGeom()
461     deleteDemandDrivenData(reconFaceCellCentresPtr_);
463     // Remote addressing and send-receive maps depend on the local
464     // position.  Therefore, it needs to be recalculated at mesh motion.
465     // Local zone addressing does not change with mesh motion
466     // HJ, 23/Jun/2011
467     deleteDemandDrivenData(remoteZoneAddressingPtr_);
469     deleteDemandDrivenData(receiveAddrPtr_);
470     deleteDemandDrivenData(sendAddrPtr_);
474 void Foam::ggiPolyPatch::clearOut()
476     clearGeom();
478     shadowIndex_ = -1;
479     zoneIndex_ = -1;
481     deleteDemandDrivenData(zoneAddressingPtr_);
482     deleteDemandDrivenData(patchToPatchPtr_);
483     deleteDemandDrivenData(localParallelPtr_);
487 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
489 Foam::ggiPolyPatch::ggiPolyPatch
491     const word& name,
492     const label size,
493     const label start,
494     const label index,
495     const polyBoundaryMesh& bm
498     coupledPolyPatch(name, size, start, index, bm),
499     shadowName_("initializeMe"),
500     zoneName_("initializeMe"),
501     bridgeOverlap_(false),
502     reject_(ggiZoneInterpolation::BB_OCTREE),
503     shadowIndex_(-1),
504     zoneIndex_(-1),
505     patchToPatchPtr_(NULL),
506     zoneAddressingPtr_(NULL),
507     remoteZoneAddressingPtr_(NULL),
508     reconFaceCellCentresPtr_(NULL),
509     localParallelPtr_(NULL),
510     receiveAddrPtr_(NULL),
511     sendAddrPtr_(NULL)
515 Foam::ggiPolyPatch::ggiPolyPatch
517     const word& name,
518     const label size,
519     const label start,
520     const label index,
521     const polyBoundaryMesh& bm,
522     const word& shadowName,
523     const word& zoneName,
524     const bool bridgeOverlap,
525     const ggiZoneInterpolation::quickReject reject
528     coupledPolyPatch(name, size, start, index, bm),
529     shadowName_(shadowName),
530     zoneName_(zoneName),
531     bridgeOverlap_(bridgeOverlap),
532     reject_(reject),
533     shadowIndex_(-1),
534     zoneIndex_(-1),
535     patchToPatchPtr_(NULL),
536     zoneAddressingPtr_(NULL),
537     remoteZoneAddressingPtr_(NULL),
538     reconFaceCellCentresPtr_(NULL),
539     localParallelPtr_(NULL),
540     receiveAddrPtr_(NULL),
541     sendAddrPtr_(NULL)
545 Foam::ggiPolyPatch::ggiPolyPatch
547     const word& name,
548     const dictionary& dict,
549     const label index,
550     const polyBoundaryMesh& bm
553     coupledPolyPatch(name, dict, index, bm),
554     shadowName_(dict.lookup("shadowPatch")),
555     zoneName_(dict.lookup("zone")),
556     bridgeOverlap_(dict.lookup("bridgeOverlap")),
557     reject_(ggiZoneInterpolation::BB_OCTREE),
558     shadowIndex_(-1),
559     zoneIndex_(-1),
560     patchToPatchPtr_(NULL),
561     zoneAddressingPtr_(NULL),
562     remoteZoneAddressingPtr_(NULL),
563     reconFaceCellCentresPtr_(NULL),
564     localParallelPtr_(NULL),
565     receiveAddrPtr_(NULL),
566     sendAddrPtr_(NULL)
568     if (dict.found("quickReject"))
569     {
570         reject_ = ggiZoneInterpolation::quickRejectNames_.read
571         (
572             dict.lookup("quickReject")
573         );
574     }
578 Foam::ggiPolyPatch::ggiPolyPatch
580     const ggiPolyPatch& pp,
581     const polyBoundaryMesh& bm
584     coupledPolyPatch(pp, bm),
585     shadowName_(pp.shadowName_),
586     zoneName_(pp.zoneName_),
587     bridgeOverlap_(pp.bridgeOverlap_),
588     reject_(pp.reject_),
589     shadowIndex_(-1),
590     zoneIndex_(-1),
591     patchToPatchPtr_(NULL),
592     zoneAddressingPtr_(NULL),
593     remoteZoneAddressingPtr_(NULL),
594     reconFaceCellCentresPtr_(NULL),
595     localParallelPtr_(NULL),
596     receiveAddrPtr_(NULL),
597     sendAddrPtr_(NULL)
601 //- Construct as copy, resetting the face list and boundary mesh data
602 Foam::ggiPolyPatch::ggiPolyPatch
604     const ggiPolyPatch& pp,
605     const polyBoundaryMesh& bm,
606     const label index,
607     const label newSize,
608     const label newStart
611     coupledPolyPatch(pp, bm, index, newSize, newStart),
612     shadowName_(pp.shadowName_),
613     zoneName_(pp.zoneName_),
614     bridgeOverlap_(pp.bridgeOverlap_),
615     reject_(pp.reject_),
616     shadowIndex_(-1),
617     zoneIndex_(-1),
618     patchToPatchPtr_(NULL),
619     zoneAddressingPtr_(NULL),
620     remoteZoneAddressingPtr_(NULL),
621     reconFaceCellCentresPtr_(NULL),
622     localParallelPtr_(NULL),
623     receiveAddrPtr_(NULL),
624     sendAddrPtr_(NULL)
628 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
630 Foam::ggiPolyPatch::~ggiPolyPatch()
632     clearOut();
636 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
638 Foam::label Foam::ggiPolyPatch::shadowIndex() const
640     if (shadowIndex_ == -1 && shadowName_ != Foam::word::null)
641     {
642         // Grab shadow patch index
643         polyPatchID shadow(shadowName_, boundaryMesh());
645         if (!shadow.active())
646         {
647             FatalErrorIn("label ggiPolyPatch::shadowIndex() const")
648                 << "Shadow patch name " << shadowName_
649                 << " not found.  Please check your GGI interface definition."
650                 << abort(FatalError);
651         }
653         shadowIndex_ = shadow.index();
655         // Check the other side is a ggi
656         if (!isA<ggiPolyPatch>(boundaryMesh()[shadowIndex_]))
657         {
658             FatalErrorIn("label ggiPolyPatch::shadowIndex() const")
659                 << "Shadow of ggi patch " << name()
660                 << " named " << shadowName() << " is not a ggi.  Type: "
661                 << boundaryMesh()[shadowIndex_].type() << nl
662                 << "This is not allowed.  Please check your mesh definition."
663                 << abort(FatalError);
664         }
666         // Check for GGI onto self
667         if (index() == shadowIndex_)
668         {
669             FatalErrorIn("label ggiPolyPatch::shadowIndex() const")
670                 << "ggi patch " << name() << " created as its own shadow"
671                 << abort(FatalError);
672         }
673     }
675     return shadowIndex_;
679 Foam::label Foam::ggiPolyPatch::zoneIndex() const
681     if (zoneIndex_ == -1 && zoneName_ != Foam::word::null)
682     {
683         // Grab zone patch index
684         faceZoneID zone(zoneName_, boundaryMesh().mesh().faceZones());
686         if (!zone.active())
687         {
688             FatalErrorIn("label ggiPolyPatch::zoneIndex() const")
689                 << "Face zone name " << zoneName_
690                 << " for GGI patch " << name()
691                 << " not found.  Please check your GGI interface definition."
692                 << abort(FatalError);
693         }
695         zoneIndex_ = zone.index();
696     }
698     return zoneIndex_;
702 const Foam::ggiPolyPatch& Foam::ggiPolyPatch::shadow() const
704     return refCast<const ggiPolyPatch>(boundaryMesh()[shadowIndex()]);
708 const Foam::faceZone& Foam::ggiPolyPatch::zone() const
710     return boundaryMesh().mesh().faceZones()[zoneIndex()];
714 const Foam::labelList& Foam::ggiPolyPatch::zoneAddressing() const
716     if (!zoneAddressingPtr_)
717     {
718         calcZoneAddressing();
719     }
721     return *zoneAddressingPtr_;
725 const Foam::labelList& Foam::ggiPolyPatch::remoteZoneAddressing() const
727     if (!remoteZoneAddressingPtr_)
728     {
729         calcRemoteZoneAddressing();
730     }
732     return *remoteZoneAddressingPtr_;
736 bool Foam::ggiPolyPatch::localParallel() const
738     // Calculate patch-to-zone addressing
739     if (!localParallelPtr_)
740     {
741         calcLocalParallel();
742     }
744     return *localParallelPtr_;
748 const Foam::ggiZoneInterpolation& Foam::ggiPolyPatch::patchToPatch() const
750     if (master())
751     {
752         if (!patchToPatchPtr_)
753         {
754             Info<< "Initializing the GGI interpolator between "
755                 << "master/shadow patches: "
756                 << name() << "/" << shadowName()
757                 << endl;
759             calcPatchToPatch();
760         }
762         return *patchToPatchPtr_;
763     }
764     else
765     {
766         return shadow().patchToPatch();
767     }
771 const Foam::vectorField& Foam::ggiPolyPatch::reconFaceCellCentres() const
773     if (!reconFaceCellCentresPtr_)
774     {
775         calcReconFaceCellCentres();
776     }
778     return *reconFaceCellCentresPtr_;
782 void Foam::ggiPolyPatch::initAddressing()
784     if (active())
785     {
786         // Calculate transforms for correct GGI cut
787         calcTransforms();
789         if (master())
790         {
791             shadow().calcTransforms();
792         }
794         // Force zone addressing and remote zone addressing
795         // (uses GGI interpolator)
796         zoneAddressing();
797         remoteZoneAddressing();
799         // Force local parallel
800         if (Pstream::parRun() && !localParallel())
801         {
802             // Calculate send addressing
803             sendAddr();
804         }
805     }
807     polyPatch::initAddressing();
811 void Foam::ggiPolyPatch::calcAddressing()
813     polyPatch::calcAddressing();
817 void Foam::ggiPolyPatch::initGeometry()
819     // Communication is allowed either before or after processor
820     // patch comms.  HJ, 11/Jul/2011
821     if (active())
822     {
823         // Note: Only master calculates recon; slave uses master interpolation
824         if (master())
825         {
826             reconFaceCellCentres();
827         }
828     }
830     polyPatch::initGeometry();
834 void Foam::ggiPolyPatch::calcGeometry()
836     polyPatch::calcGeometry();
838     // Note: Calculation of transforms must be forced before the
839     // reconFaceCellCentres in order to correctly set the transformation
840     // in the interpolation routines.
841     // HJ, 3/Jul/2009
845 void Foam::ggiPolyPatch::initMovePoints(const pointField& p)
847     clearGeom();
849     // Calculate transforms on mesh motion?
850     calcTransforms();
852     // Update interpolation for new relative position of GGI interfaces
853     if (patchToPatchPtr_)
854     {
855         patchToPatchPtr_->movePoints
856         (
857             forwardT(),
858             reverseT(),
859             -separation()
860         );
861     }
863     // Recalculate send and receive maps
864     if (active())
865     {
866         // Force zone addressing first
867         zoneAddressing();
868         remoteZoneAddressing();
870         if (Pstream::parRun() && !localParallel())
871         {
872             sendAddr();
873         }
874     }
876     if (active() && master())
877     {
878         reconFaceCellCentres();
879     }
881     polyPatch::initMovePoints(p);
885 void Foam::ggiPolyPatch::movePoints(const pointField& p)
887     polyPatch::movePoints(p);
891 void Foam::ggiPolyPatch::initUpdateMesh()
893     polyPatch::initUpdateMesh();
897 void Foam::ggiPolyPatch::updateMesh()
899     polyPatch::updateMesh();
900     clearOut();
904 void Foam::ggiPolyPatch::calcTransforms() const
906     // Simplest mixing plane: no transform or separation.  HJ, 24/Oct/2008
907     forwardT_.setSize(0);
908     reverseT_.setSize(0);
909     separation_.setSize(0);
911     if (debug > 1 && master())
912     {
913         if
914         (
915             !empty()
916          && patchToPatch().uncoveredMasterFaces().size() > 0
917         )
918         {
919             // Write uncovered master faces
920             Info<< "Writing uncovered master faces for patch "
921                 << name() << " as VTK." << endl;
923             const polyMesh& mesh = boundaryMesh().mesh();
925             fileName fvPath(mesh.time().path()/"VTK");
926             mkDir(fvPath);
928             indirectPrimitivePatch::writeVTK
929             (
930                 fvPath/fileName("uncoveredGgiFaces" + name()),
931                 IndirectList<face>
932                 (
933                     localFaces(),
934                     patchToPatch().uncoveredMasterFaces()
935                 ),
936                 localPoints()
937             );
938         }
940         if
941         (
942             !shadow().empty()
943          && patchToPatch().uncoveredSlaveFaces().size() > 0
944         )
945         {
946             // Write uncovered master faces
947             Info<< "Writing uncovered shadow faces for patch "
948                 << shadowName() << " as VTK." << endl;
950             const polyMesh& mesh = boundaryMesh().mesh();
952             fileName fvPath(mesh.time().path()/"VTK");
953             mkDir(fvPath);
954             Pout<< "shadow().localFaces(): " << shadow().localFaces().size()
955                 << " patchToPatch().uncoveredSlaveFaces().size(): "
956                 << patchToPatch().uncoveredSlaveFaces().size()
957                 << " shadow().localPoints(): " << shadow().localPoints().size()
958                 << endl;
959             indirectPrimitivePatch::writeVTK
960             (
961                 fvPath/fileName("uncoveredGgiFaces" + shadowName()),
962                 IndirectList<face>
963                 (
964                     shadow().localFaces(),
965                     patchToPatch().uncoveredSlaveFaces()
966                 ),
967                 shadow().localPoints()
968             );
969         }
971         // Check for bridge overlap
972         if (!bridgeOverlap())
973         {
974             if
975             (
976                 (
977                     patchToPatch().uncoveredMasterFaces().size() > 0
978                     && !empty()
979                 )
980              || (
981                     !shadow().empty()
982                  && patchToPatch().uncoveredSlaveFaces().size() > 0
983                 )
984             )
985             {
986                 FatalErrorIn("label ggiPolyPatch::calcTransforms() const")
987                     << "ggi patch " << name() << " with shadow "
988                     << shadowName() << " has "
989                     << patchToPatch().uncoveredMasterFaces().size()
990                     << " uncovered master faces and "
991                     << patchToPatch().uncoveredSlaveFaces().size()
992                     << " uncovered slave faces.  Bridging is switched off. "
993                     << abort(FatalError);
994             }
995         }
996         else
997         {
998             InfoIn("label ggiPolyPatch::calcTransforms() const")
999                 << "ggi patch " << name() << " with shadow "
1000                 << shadowName() << " has "
1001                 << patchToPatch().uncoveredMasterFaces().size()
1002                 << " uncovered master faces and "
1003                 << patchToPatch().uncoveredSlaveFaces().size()
1004                 << " uncovered slave faces.  Bridging is switched on. "
1005                 << endl;
1006         }
1007     }
1011 void Foam::ggiPolyPatch::initOrder(const primitivePatch&) const
1015 bool Foam::ggiPolyPatch::order
1017     const primitivePatch& pp,
1018     labelList& faceMap,
1019     labelList& rotation
1020 ) const
1022     faceMap.setSize(pp.size(), -1);
1023     rotation.setSize(pp.size(), 0);
1025     // Nothing changes
1026     return false;
1030 void Foam::ggiPolyPatch::syncOrder() const
1034 void Foam::ggiPolyPatch::write(Ostream& os) const
1036     polyPatch::write(os);
1037     os.writeKeyword("shadowPatch") << shadowName_
1038         << token::END_STATEMENT << nl;
1039     os.writeKeyword("zone") << zoneName_
1040         << token::END_STATEMENT << nl;
1041     os.writeKeyword("bridgeOverlap") << bridgeOverlap_
1042         << token::END_STATEMENT << nl;
1046 // ************************************************************************* //