Merge branch 'HrvojeJasak' of ssh://git.code.sf.net/u/hjasak/foam-extend-3.1 into...
[foam-extend-3.2.git] / src / foam / meshes / polyMesh / polyPatches / constraint / ggi / ggiPolyPatch.C
blob23deecb82073b5484d3da78bb6b586c354e6d4bb
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" << endl;
199         }
201         // Create interpolation for zones
202         patchToPatchPtr_ =
203             new ggiZoneInterpolation
204             (
205                 zone()(),           // This zone reference
206                 shadow().zone()(),  // Shadow zone reference
207                 forwardT(),
208                 reverseT(),
209                 -separation(), // Slave-to-master separation: Use - localValue
210                 // Bug fix, delayed slave evaluation causes error
211                 // HJ, 30/Jun/2013
212                 0,             // Non-overlapping face tolerances
213                 0,             // HJ, 24/Oct/2008
214                 true,          // Rescale weighting factors.  Bug fix, MB.
215                 reject_        // Quick rejection algorithm, default BB_OCTREE
216             );
218         // Abort immediately if uncovered faces are present and the option
219         // bridgeOverlap is not set.
220         if
221         (
222             (
223                 patchToPatch().uncoveredMasterFaces().size() > 0
224             && !bridgeOverlap()
225             )
226          || (
227                 patchToPatch().uncoveredSlaveFaces().size() > 0
228             && !shadow().bridgeOverlap()
229             )
230         )
231         {
232             FatalErrorIn("void ggiPolyPatch::calcPatchToPatch() const")
233                 << "Found uncovered faces for GGI interface "
234                 << name() << "/" << shadowName()
235                 << " while the bridgeOverlap option is not set "
236                 << "in the boundary file." << endl
237                 << "This is an unrecoverable error. Aborting."
238                 << abort(FatalError);
239         }
240     }
241     else
242     {
243         FatalErrorIn("void ggiPolyPatch::calcPatchToPatch() const")
244             << "Attempting to create GGIInterpolation on a shadow"
245             << abort(FatalError);
246     }
250 void Foam::ggiPolyPatch::calcReconFaceCellCentres() const
252     if (reconFaceCellCentresPtr_)
253     {
254         FatalErrorIn
255         (
256             "void ggiPolyPatch::calcReconFaceCellCentres() const"
257         )   << "Reconstructed cell centres already calculated"
258             << abort(FatalError);
259     }
261     // Create neighbouring face centres using interpolation
262     if (master())
263     {
264         const label shadowID = shadowIndex();
266         // Get the transformed and interpolated shadow face cell centers
267         reconFaceCellCentresPtr_ =
268             new vectorField
269             (
270                 interpolate
271                 (
272                     boundaryMesh()[shadowID].faceCellCentres()
273                   - boundaryMesh()[shadowID].faceCentres()
274                 )
275               + faceCentres()
276             );
277     }
278     else
279     {
280         FatalErrorIn("void ggiPolyPatch::calcReconFaceCellCentres() const")
281             << "Attempting to create reconFaceCellCentres on a shadow"
282             << abort(FatalError);
283     }
287 void Foam::ggiPolyPatch::calcLocalParallel() const
289     // Calculate patch-to-zone addressing
290     if (localParallelPtr_)
291     {
292         FatalErrorIn("void ggiPolyPatch::calcLocalParallel() const")
293             << "Local parallel switch already calculated"
294             << abort(FatalError);
295     }
297     localParallelPtr_ = new bool(false);
298     bool& emptyOrComplete = *localParallelPtr_;
300     // If running in serial, all GGIs are expanded to zone size.
301     // This happens on decomposition and reconstruction where
302     // size and shadow size may be zero, but zone size may not
303     // HJ, 1/Jun/2011
304     if (!Pstream::parRun())
305     {
306         emptyOrComplete = false;
307     }
308     else
309     {
310         // Check that patch size is greater than the zone size.
311         // This is an indication of the error where the face zone is not global
312         // in a parallel run.  HJ, 9/Nov/2014
313         if (size() > zone().size())
314         {
315             FatalErrorIn("void ggiPolyPatch::calcLocalParallel() const")
316                 << "Patch size is greater than zone size for GGI patch "
317                 << name() << ".  This is not allowerd: "
318                 << "the face zone must contain all patch faces and be "
319                 << "global in parallel runs"
320                 << abort(FatalError);
321         }
323         // Calculate localisation on master and shadow
324         emptyOrComplete =
325             (
326                 zone().size() == size()
327              && shadow().zone().size() == shadow().size()
328             )
329          || (size() == 0 && shadow().size() == 0);
331         reduce(emptyOrComplete, andOp<bool>());
332     }
334     if (debug && Pstream::parRun())
335     {
336         Info<< "GGI patch Master: " << name()
337             << " Slave: " << shadowName() << " is ";
339         if (emptyOrComplete)
340         {
341            Info<< "local parallel" << endl;
342         }
343         else
344         {
345             Info<< "split between multiple processors" << endl;
346         }
347     }
351 void Foam::ggiPolyPatch::calcSendReceive() const
353     // Note: all processors will execute calcSendReceive but only master will
354     // hold the information.  Therefore, pointers on slave processors
355     // will remain meaningless, but for purposes of consistency
356     // (of the calc-call) they will be set to zero-sized array
357     // HJ, 4/Jun/2011
359     if (receiveAddrPtr_ || sendAddrPtr_)
360     {
361         FatalErrorIn("void ggiPolyPatch::calcSendReceive() const")
362             << "Send-receive addressing already calculated"
363             << abort(FatalError);
364     }
366     if (debug)
367     {
368         Pout<< "ggiPolyPatch::calcSendReceive() const for patch "
369             << index() << endl;
370     }
372     if (!Pstream::parRun())
373     {
374         FatalErrorIn("void ggiPolyPatch::calcSendReceive() const")
375             << "Requested calculation of send-receive addressing for a "
376             << "serial run.  This is not allowed"
377             << abort(FatalError);
378     }
380     // Master will receive and store the maps
381     if (Pstream::master())
382     {
383         receiveAddrPtr_ = new labelListList(Pstream::nProcs());
384         labelListList& rAddr = *receiveAddrPtr_;
386         sendAddrPtr_ = new labelListList(Pstream::nProcs());
387         labelListList& sAddr = *sendAddrPtr_;
389         // Insert master
390         rAddr[0] = zoneAddressing();
392         for (label procI = 1; procI < Pstream::nProcs(); procI++)
393         {
394             // Note: must use normal comms because the size of the
395             // communicated lists is unknown on the receiving side
396             // HJ, 4/Jun/2011
398             // Opt: reconsider mode of communication
399             IPstream ip(Pstream::scheduled, procI);
401             rAddr[procI] = labelList(ip);
403             sAddr[procI] = labelList(ip);
404         }
405     }
406     else
407     {
408         // Create dummy pointers: only master processor stores maps
409         receiveAddrPtr_ = new labelListList();
410         sendAddrPtr_ = new labelListList();
412         // Send information to master
413         const labelList& za = zoneAddressing();
414         const labelList& ra = remoteZoneAddressing();
416         // Note: must use normal comms because the size of the
417         // communicated lists is unknown on the receiving side
418         // HJ, 4/Jun/2011
420         // Opt: reconsider mode of communication
421         OPstream op(Pstream::scheduled, Pstream::masterNo());
423         // Send local and remote addressing to master
424         op << za << ra;
425     }
429 const Foam::labelListList& Foam::ggiPolyPatch::receiveAddr() const
431     if (!receiveAddrPtr_)
432     {
433         calcSendReceive();
434     }
436     return *receiveAddrPtr_;
440 const Foam::labelListList& Foam::ggiPolyPatch::sendAddr() const
442     if (!sendAddrPtr_)
443     {
444         calcSendReceive();
445     }
447     return *sendAddrPtr_;
451 void Foam::ggiPolyPatch::clearGeom()
453     deleteDemandDrivenData(reconFaceCellCentresPtr_);
455     // Remote addressing and send-receive maps depend on the local
456     // position.  Therefore, it needs to be recalculated at mesh motion.
457     // Local zone addressing does not change with mesh motion
458     // HJ, 23/Jun/2011
459     deleteDemandDrivenData(remoteZoneAddressingPtr_);
461     deleteDemandDrivenData(receiveAddrPtr_);
462     deleteDemandDrivenData(sendAddrPtr_);
466 void Foam::ggiPolyPatch::clearOut()
468     clearGeom();
470     shadowIndex_ = -1;
471     zoneIndex_ = -1;
473     deleteDemandDrivenData(zoneAddressingPtr_);
474     deleteDemandDrivenData(patchToPatchPtr_);
475     deleteDemandDrivenData(localParallelPtr_);
479 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
481 Foam::ggiPolyPatch::ggiPolyPatch
483     const word& name,
484     const label size,
485     const label start,
486     const label index,
487     const polyBoundaryMesh& bm
490     coupledPolyPatch(name, size, start, index, bm),
491     shadowName_("initializeMe"),
492     zoneName_("initializeMe"),
493     bridgeOverlap_(false),
494     reject_(ggiZoneInterpolation::BB_OCTREE),
495     shadowIndex_(-1),
496     zoneIndex_(-1),
497     patchToPatchPtr_(NULL),
498     zoneAddressingPtr_(NULL),
499     remoteZoneAddressingPtr_(NULL),
500     reconFaceCellCentresPtr_(NULL),
501     localParallelPtr_(NULL),
502     receiveAddrPtr_(NULL),
503     sendAddrPtr_(NULL)
507 Foam::ggiPolyPatch::ggiPolyPatch
509     const word& name,
510     const label size,
511     const label start,
512     const label index,
513     const polyBoundaryMesh& bm,
514     const word& shadowName,
515     const word& zoneName,
516     const bool bridgeOverlap,
517     const ggiZoneInterpolation::quickReject reject
520     coupledPolyPatch(name, size, start, index, bm),
521     shadowName_(shadowName),
522     zoneName_(zoneName),
523     bridgeOverlap_(bridgeOverlap),
524     reject_(reject),
525     shadowIndex_(-1),
526     zoneIndex_(-1),
527     patchToPatchPtr_(NULL),
528     zoneAddressingPtr_(NULL),
529     remoteZoneAddressingPtr_(NULL),
530     reconFaceCellCentresPtr_(NULL),
531     localParallelPtr_(NULL),
532     receiveAddrPtr_(NULL),
533     sendAddrPtr_(NULL)
537 Foam::ggiPolyPatch::ggiPolyPatch
539     const word& name,
540     const dictionary& dict,
541     const label index,
542     const polyBoundaryMesh& bm
545     coupledPolyPatch(name, dict, index, bm),
546     shadowName_(dict.lookup("shadowPatch")),
547     zoneName_(dict.lookup("zone")),
548     bridgeOverlap_(dict.lookup("bridgeOverlap")),
549     reject_(ggiZoneInterpolation::BB_OCTREE),
550     shadowIndex_(-1),
551     zoneIndex_(-1),
552     patchToPatchPtr_(NULL),
553     zoneAddressingPtr_(NULL),
554     remoteZoneAddressingPtr_(NULL),
555     reconFaceCellCentresPtr_(NULL),
556     localParallelPtr_(NULL),
557     receiveAddrPtr_(NULL),
558     sendAddrPtr_(NULL)
560     if (dict.found("quickReject"))
561     {
562         reject_ = ggiZoneInterpolation::quickRejectNames_.read
563         (
564             dict.lookup("quickReject")
565         );
566     }
570 Foam::ggiPolyPatch::ggiPolyPatch
572     const ggiPolyPatch& pp,
573     const polyBoundaryMesh& bm
576     coupledPolyPatch(pp, bm),
577     shadowName_(pp.shadowName_),
578     zoneName_(pp.zoneName_),
579     bridgeOverlap_(pp.bridgeOverlap_),
580     reject_(pp.reject_),
581     shadowIndex_(-1),
582     zoneIndex_(-1),
583     patchToPatchPtr_(NULL),
584     zoneAddressingPtr_(NULL),
585     remoteZoneAddressingPtr_(NULL),
586     reconFaceCellCentresPtr_(NULL),
587     localParallelPtr_(NULL),
588     receiveAddrPtr_(NULL),
589     sendAddrPtr_(NULL)
593 //- Construct as copy, resetting the face list and boundary mesh data
594 Foam::ggiPolyPatch::ggiPolyPatch
596     const ggiPolyPatch& pp,
597     const polyBoundaryMesh& bm,
598     const label index,
599     const label newSize,
600     const label newStart
603     coupledPolyPatch(pp, bm, index, newSize, newStart),
604     shadowName_(pp.shadowName_),
605     zoneName_(pp.zoneName_),
606     bridgeOverlap_(pp.bridgeOverlap_),
607     reject_(pp.reject_),
608     shadowIndex_(-1),
609     zoneIndex_(-1),
610     patchToPatchPtr_(NULL),
611     zoneAddressingPtr_(NULL),
612     remoteZoneAddressingPtr_(NULL),
613     reconFaceCellCentresPtr_(NULL),
614     localParallelPtr_(NULL),
615     receiveAddrPtr_(NULL),
616     sendAddrPtr_(NULL)
620 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
622 Foam::ggiPolyPatch::~ggiPolyPatch()
624     clearOut();
628 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
630 Foam::label Foam::ggiPolyPatch::shadowIndex() const
632     if (shadowIndex_ == -1 && shadowName_ != Foam::word::null)
633     {
634         // Grab shadow patch index
635         polyPatchID shadow(shadowName_, boundaryMesh());
637         if (!shadow.active())
638         {
639             FatalErrorIn("label ggiPolyPatch::shadowIndex() const")
640                 << "Shadow patch name " << shadowName_
641                 << " not found.  Please check your GGI interface definition."
642                 << abort(FatalError);
643         }
645         shadowIndex_ = shadow.index();
647         // Check the other side is a ggi
648         if (!isA<ggiPolyPatch>(boundaryMesh()[shadowIndex_]))
649         {
650             FatalErrorIn("label ggiPolyPatch::shadowIndex() const")
651                 << "Shadow of ggi patch " << name()
652                 << " named " << shadowName() << " is not a ggi.  Type: "
653                 << boundaryMesh()[shadowIndex_].type() << nl
654                 << "This is not allowed.  Please check your mesh definition."
655                 << abort(FatalError);
656         }
658         // Check for GGI onto self
659         if (index() == shadowIndex_)
660         {
661             FatalErrorIn("label ggiPolyPatch::shadowIndex() const")
662                 << "ggi patch " << name() << " created as its own shadow"
663                 << abort(FatalError);
664         }
665     }
667     return shadowIndex_;
671 Foam::label Foam::ggiPolyPatch::zoneIndex() const
673     if (zoneIndex_ == -1 && zoneName_ != Foam::word::null)
674     {
675         // Grab zone patch index
676         faceZoneID zone(zoneName_, boundaryMesh().mesh().faceZones());
678         if (!zone.active())
679         {
680             FatalErrorIn("label ggiPolyPatch::zoneIndex() const")
681                 << "Face zone name " << zoneName_
682                 << " for GGI patch " << name()
683                 << " not found.  Please check your GGI interface definition."
684                 << abort(FatalError);
685         }
687         zoneIndex_ = zone.index();
688     }
690     return zoneIndex_;
694 const Foam::ggiPolyPatch& Foam::ggiPolyPatch::shadow() const
696     return refCast<const ggiPolyPatch>(boundaryMesh()[shadowIndex()]);
700 const Foam::faceZone& Foam::ggiPolyPatch::zone() const
702     return boundaryMesh().mesh().faceZones()[zoneIndex()];
706 const Foam::labelList& Foam::ggiPolyPatch::zoneAddressing() const
708     if (!zoneAddressingPtr_)
709     {
710         calcZoneAddressing();
711     }
713     return *zoneAddressingPtr_;
717 const Foam::labelList& Foam::ggiPolyPatch::remoteZoneAddressing() const
719     if (!remoteZoneAddressingPtr_)
720     {
721         calcRemoteZoneAddressing();
722     }
724     return *remoteZoneAddressingPtr_;
728 bool Foam::ggiPolyPatch::localParallel() const
730     // Calculate patch-to-zone addressing
731     if (!localParallelPtr_)
732     {
733         calcLocalParallel();
734     }
736     return *localParallelPtr_;
740 const Foam::ggiZoneInterpolation& Foam::ggiPolyPatch::patchToPatch() const
742     if (master())
743     {
744         if (!patchToPatchPtr_)
745         {
746             Info<< "Initializing the GGI interpolator between "
747                 << "master/shadow patches: "
748                 << name() << "/" << shadowName()
749                 << endl;
751             calcPatchToPatch();
752         }
754         return *patchToPatchPtr_;
755     }
756     else
757     {
758         return shadow().patchToPatch();
759     }
763 const Foam::vectorField& Foam::ggiPolyPatch::reconFaceCellCentres() const
765     if (!reconFaceCellCentresPtr_)
766     {
767         calcReconFaceCellCentres();
768     }
770     return *reconFaceCellCentresPtr_;
774 void Foam::ggiPolyPatch::initAddressing()
776     if (active())
777     {
778         // Calculate transforms for correct GGI cut
779         calcTransforms();
781         if (master())
782         {
783             shadow().calcTransforms();
784         }
786         // Force zone addressing and remote zone addressing
787         // (uses GGI interpolator)
788         zoneAddressing();
789         remoteZoneAddressing();
791         // Force local parallel
792         if (Pstream::parRun() && !localParallel())
793         {
794             // Calculate send addressing
795             sendAddr();
796         }
797     }
799     polyPatch::initAddressing();
803 void Foam::ggiPolyPatch::calcAddressing()
805     polyPatch::calcAddressing();
809 void Foam::ggiPolyPatch::initGeometry()
811     // Communication is allowed either before or after processor
812     // patch comms.  HJ, 11/Jul/2011
813     if (active())
814     {
815         // Note: Only master calculates recon; slave uses master interpolation
816         if (master())
817         {
818             reconFaceCellCentres();
819         }
820     }
822     polyPatch::initGeometry();
826 void Foam::ggiPolyPatch::calcGeometry()
828     polyPatch::calcGeometry();
830     // Note: Calculation of transforms must be forced before the
831     // reconFaceCellCentres in order to correctly set the transformation
832     // in the interpolation routines.
833     // HJ, 3/Jul/2009
837 void Foam::ggiPolyPatch::initMovePoints(const pointField& p)
839     clearGeom();
841     // Calculate transforms on mesh motion?
842     calcTransforms();
844     // Update interpolation for new relative position of GGI interfaces
845     if (patchToPatchPtr_)
846     {
847         patchToPatchPtr_->movePoints
848         (
849             forwardT(),
850             reverseT(),
851             -separation()
852         );
853     }
855     // Recalculate send and receive maps
856     if (active())
857     {
858         // Force zone addressing first
859         zoneAddressing();
860         remoteZoneAddressing();
862         if (Pstream::parRun() && !localParallel())
863         {
864             sendAddr();
865         }
866     }
868     if (active() && master())
869     {
870         reconFaceCellCentres();
871     }
873     polyPatch::initMovePoints(p);
877 void Foam::ggiPolyPatch::movePoints(const pointField& p)
879     polyPatch::movePoints(p);
883 void Foam::ggiPolyPatch::initUpdateMesh()
885     polyPatch::initUpdateMesh();
889 void Foam::ggiPolyPatch::updateMesh()
891     polyPatch::updateMesh();
892     clearOut();
896 void Foam::ggiPolyPatch::calcTransforms() const
898     // Simplest mixing plane: no transform or separation.  HJ, 24/Oct/2008
899     forwardT_.setSize(0);
900     reverseT_.setSize(0);
901     separation_.setSize(0);
903     if (debug > 1 && master())
904     {
905         if
906         (
907             !empty()
908          && patchToPatch().uncoveredMasterFaces().size() > 0
909         )
910         {
911             // Write uncovered master faces
912             Info<< "Writing uncovered master faces for patch "
913                 << name() << " as VTK." << endl;
915             const polyMesh& mesh = boundaryMesh().mesh();
917             fileName fvPath(mesh.time().path()/"VTK");
918             mkDir(fvPath);
920             indirectPrimitivePatch::writeVTK
921             (
922                 fvPath/fileName("uncoveredGgiFaces" + name()),
923                 IndirectList<face>
924                 (
925                     localFaces(),
926                     patchToPatch().uncoveredMasterFaces()
927                 ),
928                 localPoints()
929             );
930         }
932         if
933         (
934             !shadow().empty()
935          && patchToPatch().uncoveredSlaveFaces().size() > 0
936         )
937         {
938             // Write uncovered master faces
939             Info<< "Writing uncovered shadow faces for patch "
940                 << shadowName() << " as VTK." << endl;
942             const polyMesh& mesh = boundaryMesh().mesh();
944             fileName fvPath(mesh.time().path()/"VTK");
945             mkDir(fvPath);
946             Pout<< "shadow().localFaces(): " << shadow().localFaces().size()
947                 << " patchToPatch().uncoveredSlaveFaces().size(): "
948                 << patchToPatch().uncoveredSlaveFaces().size()
949                 << " shadow().localPoints(): " << shadow().localPoints().size()
950                 << endl;
951             indirectPrimitivePatch::writeVTK
952             (
953                 fvPath/fileName("uncoveredGgiFaces" + shadowName()),
954                 IndirectList<face>
955                 (
956                     shadow().localFaces(),
957                     patchToPatch().uncoveredSlaveFaces()
958                 ),
959                 shadow().localPoints()
960             );
961         }
963         // Check for bridge overlap
964         if (!bridgeOverlap())
965         {
966             if
967             (
968                 (
969                     patchToPatch().uncoveredMasterFaces().size() > 0
970                     && !empty()
971                 )
972              || (
973                     !shadow().empty()
974                  && patchToPatch().uncoveredSlaveFaces().size() > 0
975                 )
976             )
977             {
978                 FatalErrorIn("label ggiPolyPatch::calcTransforms() const")
979                     << "ggi patch " << name() << " with shadow "
980                     << shadowName() << " has "
981                     << patchToPatch().uncoveredMasterFaces().size()
982                     << " uncovered master faces and "
983                     << patchToPatch().uncoveredSlaveFaces().size()
984                     << " uncovered slave faces.  Bridging is switched off. "
985                     << abort(FatalError);
986             }
987         }
988         else
989         {
990             InfoIn("label ggiPolyPatch::calcTransforms() const")
991                 << "ggi patch " << name() << " with shadow "
992                 << shadowName() << " has "
993                 << patchToPatch().uncoveredMasterFaces().size()
994                 << " uncovered master faces and "
995                 << patchToPatch().uncoveredSlaveFaces().size()
996                 << " uncovered slave faces.  Bridging is switched on. "
997                 << endl;
998         }
999     }
1003 void Foam::ggiPolyPatch::initOrder(const primitivePatch&) const
1007 bool Foam::ggiPolyPatch::order
1009     const primitivePatch& pp,
1010     labelList& faceMap,
1011     labelList& rotation
1012 ) const
1014     faceMap.setSize(pp.size(), -1);
1015     rotation.setSize(pp.size(), 0);
1017     // Nothing changes
1018     return false;
1022 void Foam::ggiPolyPatch::syncOrder() const
1026 void Foam::ggiPolyPatch::write(Ostream& os) const
1028     polyPatch::write(os);
1029     os.writeKeyword("shadowPatch") << shadowName_
1030         << token::END_STATEMENT << nl;
1031     os.writeKeyword("zone") << zoneName_
1032         << token::END_STATEMENT << nl;
1033     os.writeKeyword("bridgeOverlap") << bridgeOverlap_
1034         << token::END_STATEMENT << nl;
1038 // ************************************************************************* //