initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / dynamicMesh / meshCut / meshModifiers / boundaryCutter / boundaryCutter.C
blobc0cfa05c13f8e8b6b67e33d4d6d0a123dc216e77
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM 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 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Description
27 \*---------------------------------------------------------------------------*/
29 #include "boundaryCutter.H"
30 #include "polyMesh.H"
31 #include "polyTopoChange.H"
32 #include "polyAddCell.H"
33 #include "polyAddFace.H"
34 #include "polyAddPoint.H"
35 #include "polyModifyFace.H"
36 #include "polyModifyPoint.H"
37 #include "mapPolyMesh.H"
38 #include "meshTools.H"
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 namespace Foam
45 defineTypeNameAndDebug(boundaryCutter, 0);
50 // * * * * * * * * * * * * * Private Static Functions  * * * * * * * * * * * //
53 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
55 void Foam::boundaryCutter::getFaceInfo
57     const label faceI,
58     label& patchID,
59     label& zoneID,
60     label& zoneFlip
61 ) const
63     patchID = -1;
65     if (!mesh_.isInternalFace(faceI))
66     {
67         patchID = mesh_.boundaryMesh().whichPatch(faceI);
68     }
70     zoneID = mesh_.faceZones().whichZone(faceI);
72     zoneFlip = false;
74     if (zoneID >= 0)
75     {
76         const faceZone& fZone = mesh_.faceZones()[zoneID];
78         zoneFlip = fZone.flipMap()[fZone.whichFace(faceI)];
79     }
83 // Adds additional vertices (from edge cutting) to face. Used for faces which
84 // are not split but still might use edge that has been cut.
85 Foam::face Foam::boundaryCutter::addEdgeCutsToFace
87     const label faceI,
88     const Map<labelList>& edgeToAddedPoints
89 ) const
91     const edgeList& edges = mesh_.edges();
92     const face& f = mesh_.faces()[faceI];
93     const labelList& fEdges = mesh_.faceEdges()[faceI];
95     // Storage for face
96     DynamicList<label> newFace(2 * f.size());
98     forAll(f, fp)
99     {
100         // Duplicate face vertex .
101         newFace.append(f[fp]);
103         // Check if edge has been cut.
104         label v1 = f.nextLabel(fp);
106         label edgeI = meshTools::findEdge(edges, fEdges, f[fp], v1);
108         Map<labelList>::const_iterator fnd = edgeToAddedPoints.find(edgeI);
110         if (fnd != edgeToAddedPoints.end())
111         {
112             // edge has been cut. Introduce new vertices. Check order.
113             const labelList& addedPoints = fnd();
115             if (edges[edgeI].start() == f[fp])
116             {
117                 // Introduce in same order.
118                 forAll(addedPoints, i)
119                 {
120                     newFace.append(addedPoints[i]);
121                 }
122             }
123             else
124             {
125                 // Introduce in opposite order.
126                 forAllReverse(addedPoints, i)
127                 {
128                     newFace.append(addedPoints[i]);
129                 }
130             }
131         }
132     }
134     face returnFace;
135     returnFace.transfer(newFace);
137     if (debug)
138     {
139         Pout<< "addEdgeCutsToFace:" << nl
140             << "    from : " << f << nl
141             << "    to   : " << returnFace << endl;
142     }
144     return returnFace;
148 void Foam::boundaryCutter::addFace
150     const label faceI,
151     const face& newFace,
153     bool& modifiedFace,     // have we already 'used' faceI
154     polyTopoChange& meshMod
155 ) const
157     // Information about old face
158     label patchID, zoneID, zoneFlip;
159     getFaceInfo(faceI, patchID, zoneID, zoneFlip);
160     label own = mesh_.faceOwner()[faceI];
161     label masterPoint = mesh_.faces()[faceI][0];
162     
163     if (!modifiedFace)
164     {
165         meshMod.setAction
166         (
167             polyModifyFace
168             (
169                 newFace,       // face
170                 faceI,
171                 own,           // owner
172                 -1,            // neighbour
173                 false,         // flux flip
174                 patchID,       // patch for face
175                 false,         // remove from zone
176                 zoneID,        // zone for face
177                 zoneFlip       // face zone flip
178             )
179         );
181         modifiedFace = true;
182     }
183     else
184     {
185         meshMod.setAction
186         (
187             polyAddFace
188             (
189                 newFace,       // face
190                 own,           // owner
191                 -1,            // neighbour
192                 masterPoint,   // master point
193                 -1,            // master edge
194                 -1,            // master face for addition
195                 false,         // flux flip
196                 patchID,       // patch for face
197                 zoneID,        // zone for face
198                 zoneFlip       // face zone flip
199             )
200         );
201     }
206 // Splits a face using the cut edges and modified points 
207 bool Foam::boundaryCutter::splitFace
209     const label faceI,
210     const Map<point>& pointToPos,
211     const Map<labelList>& edgeToAddedPoints,
212     polyTopoChange& meshMod
213 ) const
215     const edgeList& edges = mesh_.edges();
216     const face& f = mesh_.faces()[faceI];
217     const labelList& fEdges = mesh_.faceEdges()[faceI];
219     // Count number of split edges and total number of splits.
220     label nSplitEdges = 0;
221     label nModPoints = 0;
222     label nTotalSplits = 0;
224     forAll(f, fp)
225     {
226         if (pointToPos.found(f[fp]))
227         {
228             nModPoints++;
229             nTotalSplits++;
230         }
232         // Check if edge has been cut.
233         label nextV = f.nextLabel(fp);
235         label edgeI = meshTools::findEdge(edges, fEdges, f[fp], nextV);
237         Map<labelList>::const_iterator fnd = edgeToAddedPoints.find(edgeI);
239         if (fnd != edgeToAddedPoints.end())
240         {
241             nSplitEdges++;
242             nTotalSplits += fnd().size();
243         }
244     }
246     if (debug)
247     {
248         Pout<< "Face:" << faceI
249             << " nModPoints:" << nModPoints
250             << " nSplitEdges:" << nSplitEdges
251             << " nTotalSplits:" << nTotalSplits << endl;
252     }
254     if (nSplitEdges == 0 && nModPoints == 0)
255     {
256         FatalErrorIn("boundaryCutter::splitFace") << "Problem : face:" << faceI
257             << " nSplitEdges:" << nSplitEdges
258             << " nTotalSplits:" << nTotalSplits
259             << abort(FatalError);
260         return false;
261     }
262     else if (nSplitEdges + nModPoints == 1)
263     {
264         // single or multiple cuts on a single edge or single modified point
265         // Dont cut and let caller handle this.
266         Warning << "Face " << faceI << " has only one edge cut " << endl;
267         return false;
268     }
269     else
270     {
271         // So guaranteed to have two edges cut or points modified. Split face:
272         // - find starting cut
273         // - walk to next cut. Make face
274         // - loop until face done.
276         // Information about old face
277         label patchID, zoneID, zoneFlip;
278         getFaceInfo(faceI, patchID, zoneID, zoneFlip);
280         // Get face with new points on cut edges for ease of looping
281         face extendedFace(addEdgeCutsToFace(faceI, edgeToAddedPoints));
283         // Find first added point. This is the starting vertex for splitting.
284         label startFp = -1;
286         forAll(extendedFace, fp)
287         {
288             if (extendedFace[fp] >= mesh_.nPoints())
289             {
290                 startFp = fp;
291                 break;
292             }
293         }
295         if (startFp == -1)
296         {
297             // No added point. Maybe there is a modified point?
298             forAll(extendedFace, fp)
299             {
300                 if (pointToPos.found(extendedFace[fp]))
301                 {
302                     startFp = fp;
303                     break;
304                 }
305             }
306         }
308         if (startFp == -1)
309         {
310             FatalErrorIn("boundaryCutter::splitFace")
311                 << "Problem" << abort(FatalError);
312         }
314         // Have we already modified existing face (first face gets done
315         // as modification; all following ones as polyAddFace)
316         bool modifiedFace = false;
318         // Example face:
319         //    +--+
320         //   /   |
321         //  /    |
322         // +     +
323         //  \    |
324         //   \   |
325         //    +--+
326         //
327         // Needs to get split into:
328         // - three 'side' faces a,b,c
329         // - one middle face d
330         //    +--+
331         //   /|\A|
332         //  / | \|
333         // + C|D +
334         //  \ | /|
335         //   \|/B|
336         //    +--+
339         // Storage for new face
340         DynamicList<label> newFace(extendedFace.size());
342         label fp = startFp;
344         forAll(extendedFace, i)
345         {
346             label pointI = extendedFace[fp];
348             newFace.append(pointI);
350             if
351             (
352                 newFace.size() > 2
353              && (
354                     pointI >= mesh_.nPoints()
355                  || pointToPos.found(pointI)
356                 )
357             )
358             {
359                 // Enough vertices to create a face from.
360                 face tmpFace;
361                 tmpFace.transfer(newFace);
363                 // Add face tmpFace
364                 addFace(faceI, tmpFace, modifiedFace, meshMod);
366                 // Starting point is also the starting point for the new face
367                 newFace.append(extendedFace[startFp]);
368                 newFace.append(extendedFace[fp]);
369             }
371             fp = (fp+1) % extendedFace.size();
372         }
374         // Check final face.
375         if (newFace.size() > 2)
376         {
377             // Enough vertices to create a face from.
378             face tmpFace;
379             tmpFace.transfer(newFace);
381             // Add face tmpFace
382             addFace(faceI, tmpFace, modifiedFace, meshMod);
383         }
385         // Split something
386         return true;
387     }
391 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
393 // Construct from components
394 Foam::boundaryCutter::boundaryCutter(const polyMesh& mesh)
396     mesh_(mesh),
397     edgeAddedPoints_(),
398     faceAddedPoint_()
402 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
404 Foam::boundaryCutter::~boundaryCutter()
408 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
410 void Foam::boundaryCutter::setRefinement
412     const Map<point>& pointToPos,
413     const Map<List<point> >& edgeToCuts,
414     const Map<labelPair>& faceToSplit,
415     const Map<point>& faceToFeaturePoint,
416     polyTopoChange& meshMod
419     // Clear and size maps here since mesh size will change.
420     edgeAddedPoints_.clear();
422     faceAddedPoint_.clear();
423     faceAddedPoint_.resize(faceToFeaturePoint.size());
426     //
427     // Points that just need to be moved
428     // Note: could just as well be handled outside of setRefinement.
429     //
431     forAllConstIter(Map<point>, pointToPos, iter)
432     {
433         meshMod.setAction
434         (
435             polyModifyPoint
436             (
437                 iter.key(), // point
438                 iter(),     // position
439                 false,      // no zone
440                 -1,         // zone for point
441                 true        // supports a cell
442             )
443         );
444     }
447     //
448     // Add new points along cut edges.
449     //
451     // Map from edge label to sorted list of points
452     Map<labelList> edgeToAddedPoints(edgeToCuts.size());
454     forAllConstIter(Map<List<point> >, edgeToCuts, iter)
455     {
456         label edgeI = iter.key();
458         const edge& e = mesh_.edges()[edgeI];
460         // Sorted (from start to end) list of cuts on edge
461         const List<point>& cuts = iter();
463         forAll(cuts, cutI)
464         {
465             // point on feature to move to
466             const point& featurePoint = cuts[cutI];
468             label addedPointI =
469                 meshMod.setAction
470                 (
471                     polyAddPoint
472                     (
473                         featurePoint,               // point
474                         e.start(),                  // master point
475                         -1,                         // zone for point
476                         true                        // supports a cell
477                     )
478                 );
480             Map<labelList>::iterator fnd = edgeToAddedPoints.find(edgeI);
482             if (fnd != edgeToAddedPoints.end())
483             {
484                 labelList& addedPoints = fnd();
486                 label sz = addedPoints.size();
487                 addedPoints.setSize(sz+1);
488                 addedPoints[sz] = addedPointI;
489             }
490             else
491             {
492                 edgeToAddedPoints.insert(edgeI, labelList(1, addedPointI));
493             }
495             if (debug)
496             {
497                 Pout<< "Added point " << addedPointI << " for edge " << edgeI
498                     << " with cuts:" << edgeToAddedPoints[edgeI] << endl;
499             }
500         }
501     }
504     //
505     // Introduce feature points.
506     //
508     forAllConstIter(Map<point>, faceToFeaturePoint, iter)
509     {
510         label faceI = iter.key();
512         const face& f = mesh_.faces()[faceI];
514         if (faceToSplit.found(faceI))
515         {
516             FatalErrorIn("boundaryCutter::setRefinement")
517                 << "Face " << faceI << " vertices " << f
518                 << " is both marked for face-centre decomposition and"
519                 << " diagonal splitting."
520                 << abort(FatalError);
521         }
523         if (mesh_.isInternalFace(faceI))
524         {
525             FatalErrorIn("boundaryCutter::setRefinement")
526                 << "Face " << faceI << " vertices " << f
527                 << " is not an external face. Cannot split it"
528                 << abort(FatalError);
529         }
531         label addedPointI =
532             meshMod.setAction
533             (
534                 polyAddPoint
535                 (
536                     iter(), // point
537                     f[0],   // master point
538                     -1,     // zone for point
539                     true    // supports a cell
540                 )
541             );
542         faceAddedPoint_.insert(faceI, addedPointI);
544         if (debug)
545         {
546             Pout<< "Added point " << addedPointI << " for feature point "
547                 << iter() << " on face " << faceI << " with centre "
548                 << mesh_.faceCentres()[faceI] << endl;
549         }
550     }
553     //
554     // Split or retriangulate faces
555     //
558     // Maintain whether face has been updated (for -split edges
559     // -new owner/neighbour)
560     boolList faceUptodate(mesh_.nFaces(), false);
563     // Triangulate faces containing feature points
564     forAllConstIter(Map<label>, faceAddedPoint_, iter)
565     {
566         label faceI = iter.key();
568         // Get face with new points on cut edges.
569         face newFace(addEdgeCutsToFace(faceI, edgeToAddedPoints));
571         label addedPointI = iter();
573         // Information about old face
574         label patchID, zoneID, zoneFlip;
575         getFaceInfo(faceI, patchID, zoneID, zoneFlip);
576         label own = mesh_.faceOwner()[faceI];
577         label masterPoint = mesh_.faces()[faceI][0];
579         // Triangulate face around mid point
581         face tri(3);
583         forAll(newFace, fp)
584         {
585             label nextV = newFace.nextLabel(fp);
587             tri[0] = newFace[fp];
588             tri[1] = nextV;
589             tri[2] = addedPointI;
591             if (fp == 0)
592             {
593                 // Modify the existing face.
594                 meshMod.setAction
595                 (
596                     polyModifyFace
597                     (
598                         tri,                        // face
599                         faceI,
600                         own,                        // owner
601                         -1,                         // neighbour
602                         false,                      // flux flip
603                         patchID,                    // patch for face
604                         false,                      // remove from zone
605                         zoneID,                     // zone for face
606                         zoneFlip                    // face zone flip
607                     )
608                 );
609             }
610             else
611             {
612                 // Add additional faces
613                 meshMod.setAction
614                 (
615                     polyAddFace
616                     (
617                         tri,                        // face
618                         own,                        // owner
619                         -1,                         // neighbour
620                         masterPoint,                // master point
621                         -1,                         // master edge
622                         -1,                         // master face for addition
623                         false,                      // flux flip
624                         patchID,                    // patch for face
625                         zoneID,                     // zone for face
626                         zoneFlip                    // face zone flip
627                     )
628                 );
629             }
630         }
632         faceUptodate[faceI] = true;
633     }
636     // Diagonally split faces
637     forAllConstIter(Map<labelPair>, faceToSplit, iter)
638     {
639         label faceI = iter.key();
641         const face& f = mesh_.faces()[faceI];
643         if (faceAddedPoint_.found(faceI))
644         {
645             FatalErrorIn("boundaryCutter::setRefinement")
646                 << "Face " << faceI << " vertices " << f
647                 << " is both marked for face-centre decomposition and"
648                 << " diagonal splitting."
649                 << abort(FatalError);
650         }
653         // Get face with new points on cut edges.
654         face newFace(addEdgeCutsToFace(faceI, edgeToAddedPoints));
656         // Information about old face
657         label patchID, zoneID, zoneFlip;
658         getFaceInfo(faceI, patchID, zoneID, zoneFlip);
659         label own = mesh_.faceOwner()[faceI];
660         label masterPoint = mesh_.faces()[faceI][0];
662         // Split face from one side of diagonal to other.
663         const labelPair& diag = iter();
665         label fp0 = findIndex(newFace, f[diag[0]]);
666         label fp1 = findIndex(newFace, f[diag[1]]);
668         if (fp0 == -1 || fp1 == -1 || fp0 == fp1)
669         {
670             FatalErrorIn("boundaryCutter::setRefinement")
671                 << "Problem : Face " << faceI << " vertices " << f
672                 << " newFace:" << newFace << " diagonal:" << f[diag[0]]
673                 << ' ' << f[diag[1]]
674                 << abort(FatalError);
675         }
677         // Replace existing face by newFace from fp0 to fp1 and add new one
678         // from fp1 to fp0.
680         DynamicList<label> newVerts(newFace.size());
682         // Get vertices from fp0 to (and including) fp1
683         label fp = fp0;
685         do
686         {
687             newVerts.append(newFace[fp]);
689             fp = (fp == newFace.size()-1 ? 0 : fp+1);
691         } while (fp != fp1);
693         newVerts.append(newFace[fp1]);
696         // Modify the existing face.
697         meshMod.setAction
698         (
699             polyModifyFace
700             (
701                 face(newVerts.shrink()),    // face
702                 faceI,
703                 own,                        // owner
704                 -1,                         // neighbour
705                 false,                      // flux flip
706                 patchID,                    // patch for face
707                 false,                      // remove from zone
708                 zoneID,                     // zone for face
709                 zoneFlip                    // face zone flip
710             )
711         );
714         newVerts.clear();
716         // Get vertices from fp1 to (and including) fp0
718         do
719         {
720             newVerts.append(newFace[fp]);
722             fp = (fp == newFace.size()-1 ? 0 : fp+1);
724         } while (fp != fp0);
726         newVerts.append(newFace[fp0]);
728         // Add additional face
729         meshMod.setAction
730         (
731             polyAddFace
732             (
733                 face(newVerts.shrink()),    // face
734                 own,                        // owner
735                 -1,                         // neighbour
736                 masterPoint,                // master point
737                 -1,                         // master edge
738                 -1,                         // master face for addition
739                 false,                      // flux flip
740                 patchID,                    // patch for face
741                 zoneID,                     // zone for face
742                 zoneFlip                    // face zone flip
743             )
744         );
746         faceUptodate[faceI] = true;
747     }
750     // Split external faces without feature point but using cut edges.
751     // Does right handed walk but not really.
752     forAllConstIter(Map<labelList>, edgeToAddedPoints, iter)
753     {
754         label edgeI = iter.key();
756         const labelList& eFaces = mesh_.edgeFaces()[edgeI];
758         forAll(eFaces, i)    
759         {
760             label faceI = eFaces[i];
762             if (!faceUptodate[faceI] && !mesh_.isInternalFace(faceI))
763             {
764                 // Is external face so split
765                 if (splitFace(faceI, pointToPos, edgeToAddedPoints, meshMod))
766                 {
767                     // Successfull split
768                     faceUptodate[faceI] = true;
769                 }
770             }
771         }
772     }
775     // Add cut edges (but don't split) any other faces using any cut edge.
776     // These can be external faces where splitFace hasn't cut them or
777     // internal faces.
778     forAllConstIter(Map<labelList>, edgeToAddedPoints, iter)
779     {
780         label edgeI = iter.key();
782         const labelList& eFaces = mesh_.edgeFaces()[edgeI];
784         forAll(eFaces, i)    
785         {
786             label faceI = eFaces[i];
788             if (!faceUptodate[faceI])
789             {
790                 // Renumber face to include split edges.
791                 face newFace(addEdgeCutsToFace(faceI, edgeToAddedPoints));
793                 label own = mesh_.faceOwner()[faceI];
795                 label nei = -1;
797                 if (mesh_.isInternalFace(faceI))
798                 {
799                     nei = mesh_.faceNeighbour()[faceI];
800                 }
802                 label patchID, zoneID, zoneFlip;
803                 getFaceInfo(faceI, patchID, zoneID, zoneFlip);
805                 meshMod.setAction
806                 (
807                     polyModifyFace
808                     (
809                         newFace,            // modified face
810                         faceI,              // label of face being modified
811                         own,                // owner
812                         nei,                // neighbour
813                         false,              // face flip
814                         patchID,            // patch for face
815                         false,              // remove from zone
816                         zoneID,             // zone for face
817                         zoneFlip            // face flip in zone
818                     )
819                 );
821                 faceUptodate[faceI] = true;
822             }
823         }
824     }
826     // Convert edge to points storage from edge labels (not preserved)
827     // to point labels
828     edgeAddedPoints_.resize(edgeToCuts.size());
830     forAllConstIter(Map<labelList>, edgeToAddedPoints, iter)
831     {
832         edgeAddedPoints_.insert(mesh_.edges()[iter.key()], iter());
833     }
837 void Foam::boundaryCutter::updateMesh(const mapPolyMesh& morphMap)
839     // Update stored labels for mesh change.
841     //
842     // Do faceToAddedPoint
843     //
845     {
846         // Create copy since we're deleting entries.
847         Map<label> newAddedPoints(faceAddedPoint_.size());
849         forAllConstIter(Map<label>, faceAddedPoint_, iter)
850         {
851             label oldFaceI = iter.key();
853             label newFaceI = morphMap.reverseFaceMap()[oldFaceI];
855             label oldPointI = iter();
857             label newPointI = morphMap.reversePointMap()[oldPointI];
859             if (newFaceI >= 0 && newPointI >= 0)
860             {
861                 newAddedPoints.insert(newFaceI, newPointI);
862             }
863         }
865         // Copy
866         faceAddedPoint_.transfer(newAddedPoints);
867     }
870     //
871     // Do edgeToAddedPoints
872     //
875     {
876         // Create copy since we're deleting entries
877         HashTable<labelList, edge, Hash<edge> >
878             newEdgeAddedPoints(edgeAddedPoints_.size());
880         for
881         (
882             HashTable<labelList, edge, Hash<edge> >::const_iterator iter =
883                 edgeAddedPoints_.begin();
884             iter != edgeAddedPoints_.end();
885             ++iter
886         )
887         {
888             const edge& e = iter.key();
890             label newStart = morphMap.reversePointMap()[e.start()];
892             label newEnd = morphMap.reversePointMap()[e.end()];
894             if (newStart >= 0 && newEnd >= 0)
895             {
896                 const labelList& addedPoints = iter();
898                 labelList newAddedPoints(addedPoints.size());
899                 label newI = 0;
901                 forAll(addedPoints, i)
902                 {
903                     label newAddedPointI =
904                         morphMap.reversePointMap()[addedPoints[i]];
906                     if (newAddedPointI >= 0)
907                     {
908                         newAddedPoints[newI++] = newAddedPointI;
909                     }
910                 }
911                 if (newI > 0)
912                 {
913                     newAddedPoints.setSize(newI);
915                     edge newE = edge(newStart, newEnd);
917                     newEdgeAddedPoints.insert(newE, newAddedPoints);
918                 }
919             }
920         }
922         // Copy
923         edgeAddedPoints_.transfer(newEdgeAddedPoints);
924     }
928 // ************************************************************************* //