Fix tutorials: typo in tutorials/viscoelastic/viscoelasticFluidFoam/S-MDCPP/constant...
[OpenFOAM-1.6-ext.git] / src / dynamicMesh / polyMeshModifiers / slidingInterface / enrichedPatch / enrichedPatchCutFaces.C
blobaf510206a7660c2e41f0648cacfc06cf35dd868b
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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
26     Calculating cut faces of the enriched patch, together with the addressing
27     into master and slave patch.
29 Author
30     Hrvoje Jasak, Wikki Ltd.  All rights reserved.  Copyright Hrvoje Jasak.
32 \*---------------------------------------------------------------------------*/
34 #include "enrichedPatch.H"
35 #include "boolList.H"
36 #include "DynamicList.H"
37 #include "labelPair.H"
38 #include "primitiveMesh.H"
39 #include "HashSet.H"
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 // If the cut face gets more than this number of points, it will be checked
44 const Foam::label Foam::enrichedPatch::maxFaceSizeDebug_ = 100;
46 // Out of plane tolerace.  Used to reject points in the bent interfaces
47 const Foam::scalar Foam::enrichedPatch::outOfPlaneTol_ = 0.2;
50 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
52 // Index of debug signs:
53 // x - skip a point
54 // l - left turn
55 // r - right turn
57 void Foam::enrichedPatch::calcCutFaces() const
59     if (cutFacesPtr_ || cutFaceMasterPtr_ || cutFaceSlavePtr_)
60     {
61         FatalErrorIn("void enrichedPatch::calcCutFaces() const")
62             << "Cut faces addressing already calculated."
63             << abort(FatalError);
64     }
66     const faceList& enFaces = enrichedFaces();
67     const labelList& mp = meshPoints();
68     const faceList& lf = localFaces();
69     const pointField& lp = localPoints();
70     const labelListList& pp = pointPoints();
71 //     Pout << "enFaces: " << enFaces << endl;
72 //     Pout << "lf: " << lf << endl;
73 //     Pout << "lp: " << lp << endl;
74 //     Pout << "pp: " << pp << endl;
75     const Map<labelList>& masterPointFaceAddr = masterPointFaces();
77     // Prepare the storage
78     DynamicList<face> cf(2*lf.size());
79     DynamicList<label> cfMaster(2*lf.size());
80     DynamicList<label> cfSlave(2*lf.size());
82     // Algorithm
83     // Go through all the faces
84     // 1) For each face, start at point zero and grab the edge zero.
85     //    Both points are added into the cut face.  Mark the first edge
86     //    as used and position the current point as the end of the first
87     //    edge and previous point as point zero.
88     // 2) Grab all possible points from the current point.  Excluding
89     //    the previous point find the point giving the biggest right
90     //    turn. Add the point to the face and mark the edges as used.
91     //    Continue doing this until the face is complete, i.e. the next point
92     //    to add is the first point of the face.
93     // 3) Once the facet is completed, register its owner from the face
94     //    it has been created from (remember that the first lot of faces
95     //    inserted into the enriched faces list are the slaves, to allow
96     //    matching of the other side).
97     // 4) If the facet is created from an original slave face, go
98     //    through the master patch and try to identify the master face
99     //    this facet belongs to.  This is recognised by the fact that the
100     //    faces consists exclusively of the points of the master face and
101     //    the points projected onto the face.
103     // Create a set of edge usage parameters
104     HashSet<edge, Hash<edge> > edgesUsedOnce(pp.size());
105     HashSet<edge, Hash<edge> > edgesUsedTwice
106         (pp.size()*primitiveMesh::edgesPerPoint_);
109     forAll (lf, faceI)
110     {
111         const face& curLocalFace = lf[faceI];
112         const face& curGlobalFace = enFaces[faceI];
114 //         Pout<< "Doing face " << faceI << " local: " << curLocalFace << " or " << curGlobalFace << endl;
115 //         if (faceI < slavePatch_.size())
116 //         {
117 //             Pout<< "original slave: " << slavePatch_[faceI]
118 //                 << " local: " << slavePatch_.localFaces()[faceI] << endl;
119 //         }
120 //         else
121 //         {
122 //             Pout<< "original master: "
123 //                 << masterPatch_[faceI - slavePatch_.size()] << " "
124 //                 << masterPatch_.localFaces()[faceI - slavePatch_.size()]
125 //                 << endl;
126 //         }
127 //         {
128 //             pointField facePoints = curLocalFace.points(lp);
129 //             forAll (curLocalFace, pointI)
130 //             {
131 //                 Pout << "v " << facePoints[pointI].x() << " "
132 //                     << facePoints[pointI].y() << " "
133 //                     << facePoints[pointI].z() << endl;
134 //             }
135 //         }
137         // Track the usage of face edges.  When all edges are used, the
138         // face decomposition is complete.
139         // The seed edges include all the edges of the original face + all edges
140         // of other faces that have been used in the construction of the
141         // facet.  Edges from other faces can be considered as
142         // internal to the current face if used only once.
144         // Track the edge usage to avoid duplicate faces and reset it to unused
145         boolList usedFaceEdges(curLocalFace.size(), false);
147         SLList<edge> edgeSeeds;
149         // Insert the edges of current face into the seed list.
150         edgeList cfe = curLocalFace.edges();
151         forAll (curLocalFace, edgeI)
152         {
153             edgeSeeds.append(cfe[edgeI]);
154         }
156         // Grab face normal
157         vector normal = curLocalFace.normal(lp);
158         normal /= mag(normal);
160         while (edgeSeeds.size() > 0)
161         {
162 //             Pout << "edgeSeeds.size(): " << edgeSeeds.size() << endl;
163             const edge curEdge = edgeSeeds.removeHead();
165             // Locate the edge in current face
166             const label curEdgeWhich = curLocalFace.which(curEdge.start());
168             // Check if the edge is in current face and if it has been
169             // used already.  If so, skip it
170             if
171             (
172                 curEdgeWhich > -1
173              && curLocalFace.nextLabel(curEdgeWhich) == curEdge.end()
174             )
175             {
176                 // Edge is in the starting face.
177                 // If unused, mark it as used; if used, skip it
178                 if (usedFaceEdges[curEdgeWhich]) continue;
180                 usedFaceEdges[curEdgeWhich] = true;
181             }
183             // If the edge has already been used twice, skip it
184             if (edgesUsedTwice.found(curEdge)) continue;
185 //             Pout << "Trying new edge (" << mp[curEdge.start()] << ", " << mp[curEdge.end()] << ") seed: " << curEdge << " used: " << edgesUsedTwice.found(curEdge) << endl;
187             // Estimate the size of cut face as twice the size of original face
188             DynamicList<label> cutFaceGlobalPoints(2*curLocalFace.size());
189             DynamicList<label> cutFaceLocalPoints(2*curLocalFace.size());
191             // Found unused edge.
192             label prevPointLabel = curEdge.start();
193             cutFaceGlobalPoints.append(mp[prevPointLabel]);
194             cutFaceLocalPoints.append(prevPointLabel);
195 //              Pout << "prevPointLabel: " << mp[prevPointLabel] << endl;
196             // Grab current point and append it to the list
197             label curPointLabel = curEdge.end();
198             point curPoint = lp[curPointLabel];
199             cutFaceGlobalPoints.append(mp[curPointLabel]);
200             cutFaceLocalPoints.append(curPointLabel);
202             bool completedCutFace = false;
204             label faceSizeDebug = maxFaceSizeDebug_;
206             do
207             {
208                 // Grab the next point options
209 //                 Pout << "curPointLabel: " << mp[curPointLabel] << endl;
210                 const labelList& nextPoints = pp[curPointLabel];
211 //                 Pout << "nextPoints: " << IndirectList<label>(mp, nextPoints)() << endl;
212                 // Get the vector along the edge and the right vector
213                 vector ahead = curPoint - lp[prevPointLabel];
214                 ahead -= normal*(normal & ahead);
215                 ahead /= mag(ahead);
217                 vector right = normal ^ ahead;
218                 right /= mag(right);
219 //                 Pout<< "normal: " << normal << " ahead: " << ahead << " right: " << right << endl;
221                 scalar atanTurn = -GREAT;
222                 label bestAtanPoint = -1;
224                 forAll (nextPoints, nextI)
225                 {
226                     // Exclude the point we are coming from; there will always
227                     // be more than one edge, so this is safe
228                     if (nextPoints[nextI] != prevPointLabel)
229                     {
230 //                         Pout << "cur point: " << curPoint << " trying for point: " << mp[nextPoints[nextI]] << " " << lp[nextPoints[nextI]];
231                         vector newDir = lp[nextPoints[nextI]] - curPoint;
232                         scalar magNewDir = mag(newDir);
233 //                         Pout << " newDir: " << newDir << " mag: " << magNewDir << flush;
234                         vector corrNewDir = newDir - normal*(normal & newDir);
235                         scalar magCorrNewDir = mag(corrNewDir);
237                         // Check how far out of plane is the point
238                         scalar outOfPlane =
239                             (magNewDir - magCorrNewDir)/magNewDir;
240 //                         Pout << " corrected: " << corrNewDir << " mag: " << magCorrNewDir << " out of plane: " << outOfPlane << flush;
242                         if (magCorrNewDir < SMALL)
243                         {
244                             FatalErrorIn
245                             (
246                                 "void enrichedPatch::"
247                                 "calcCutFaces() const"
248                             )   << "Zero length edge detected.  Probable "
249                                 << "projection error: slave patch probably "
250                                 << "does not project onto master.  "
251                                 << "Please switch on "
252                                 << "enriched patch debug for more info"
253                                 << abort(FatalError);
254                         }
256                         corrNewDir /= magCorrNewDir;
258                         scalar curAtanTurn =
259                             atan2(corrNewDir & right, corrNewDir & ahead);
261 //                         Pout << " atan: " << curAtanTurn << endl;
263                         if
264                         (
265                             curAtanTurn > atanTurn
266                          && outOfPlane < outOfPlaneTol_
267                         )
268                         {
269                             bestAtanPoint = nextPoints[nextI];
270                             atanTurn = curAtanTurn;
271                         }
272                     } // end of prev point skip
273                 } // end of next point selection
275                 if (bestAtanPoint < 0)
276                 {
277                     FatalErrorIn("void enrichedPatch::calcCutFaces() const")
278                         << "Cannot determine best atan point.  "
279                         << "Doing face: " << curGlobalFace
280                         << " Cur point: " << mp[curPointLabel]
281                         << " Next points: "
282                         << IndirectList<label>(mp, nextPoints)()
283                         << abort(FatalError);
284                 }
285 //                 Pout<< "   bestAtanPoint: " << bestAtanPoint << " or "
286 //                     << mp[bestAtanPoint] << endl;
287                 // Selected next best point.
288 //                 Pout<< "cutFaceGlobalPoints: " << cutFaceGlobalPoints
289 //                     << endl;
290                 // Check if the edge about to be added has been used
291                 // in the current face or twice in other faces.  If
292                 // so, the face is bad.
293                 const label whichNextPoint = curLocalFace.which(curPointLabel);
295                 if
296                 (
297                     whichNextPoint > -1
298                  && curLocalFace.nextLabel(whichNextPoint) == bestAtanPoint
299                  && usedFaceEdges[whichNextPoint]
300                 )
301                 {
302                     // This edge is already used in current face
303                     // face cannot be good; start on a new one
304 //                     Pout<< "Double usage in current face, cannot be good"
305 //                         << endl;
306                     completedCutFace = true;
307                 }
310                 if (edgesUsedTwice.found(edge(curPointLabel, bestAtanPoint)))
311                 {
312                     // This edge is already used -
313                     // face cannot be good; start on a new one
314                     completedCutFace = true;
315 //                     Pout<< "Double usage elsewhere, cannot be good" << endl;
316                 }
318                 if (completedCutFace) continue;
320                 // Insert the next best point into the list
321                 if (mp[bestAtanPoint] == cutFaceGlobalPoints[0])
322                 {
323                     // Face is completed.  Save it and move on to the next
324                     // available edge
325                     completedCutFace = true;
327                     // Append the face
328                     face cutFaceGlobal;
329                     cutFaceGlobal.transfer(cutFaceGlobalPoints.shrink());
331                     cf.append(cutFaceGlobal);
333                     face cutFaceLocal;
334                     cutFaceLocal.transfer(cutFaceLocalPoints.shrink());
336                     if (debug)
337                     {
338                         Pout<< " local: " << cutFaceLocal
339                             << " area: " << cutFaceLocal.mag(lp)
340                             << " one side: " << faceI;
341                     }
343 //                     Pout<< "\ncutFaceLocal: " << cutFaceLocal.points(lp) << endl;
344                     // Go through all edges of the cut faces.
345                     // If the edge corresponds to a starting face edge,
346                     // mark the starting face edge as true
348                     forAll (cutFaceLocal, cutI)
349                     {
350                         const edge curCutFaceEdge
351                         (
352                             cutFaceLocal[cutI],
353                             cutFaceLocal.nextLabel(cutI)
354                         );
356                         // Increment the usage count using two hash sets
357                         HashSet<edge, Hash<edge> >::iterator euoIter =
358                             edgesUsedOnce.find(curCutFaceEdge);
360                         if (euoIter == edgesUsedOnce.end())
361                         {
362 //                             Pout<< "Found edge not used before: "
363 //                                 << curCutFaceEdge << endl;
364                             edgesUsedOnce.insert(curCutFaceEdge);
365                         }
366                         else
367                         {
368 //                             Pout<< "Found edge used once: "
369 //                                 << curCutFaceEdge << endl;
370                             edgesUsedOnce.erase(euoIter);
371                             edgesUsedTwice.insert(curCutFaceEdge);
372                         }
374                         const label curCutFaceEdgeWhich =
375                             curLocalFace.which(curCutFaceEdge.start());
377                         if
378                         (
379                             curCutFaceEdgeWhich > -1
380                          && (
381                                 curLocalFace.nextLabel(curCutFaceEdgeWhich)
382                              == curCutFaceEdge.end()
383                             )
384                         )
385                         {
386                             // Found edge in original face
387 //                             Pout<< "Found edge in orig face: "
388 //                                 << curCutFaceEdge << ": "
389 //                                 << curCutFaceEdgeWhich << endl;
390                             usedFaceEdges[curCutFaceEdgeWhich] = true;
391                         }
392                         else
393                         {
394                             // Edge not in original face.  Add it to seeds
395 //                             Pout<< "Found new edge seed: "
396 //                                 << curCutFaceEdge << endl;
397                             edgeSeeds.append(curCutFaceEdge.reverseEdge());
398                         }
399                     }
402                     // Find out what the other side is
404                     // Algorithm
405                     // If the face is in the slave half of the
406                     // enrichedFaces lists, it may be matched against
407                     // the master face.  It can be recognised by the
408                     // fact that all its points belong to one master
409                     // face.  For this purpose:
410                     // 1) Grab the master faces around the point zero
411                     // of the cut face and collect all master faces
412                     // around it.
413                     // 2) Loop through the rest of cut face points and
414                     // if the point refers to one of the faces marked
415                     // by point zero, increment its count.
416                     // 3) When completed, the face whose count is
417                     // equal to the number of points in the cut face
418                     // is the other side.  If this is not the case, there is no
419                     // face on the other side.
421                     if (faceI < slavePatch_.size())
422                     {
423                         Map<labelList>::const_iterator mpfAddrIter =
424                             masterPointFaceAddr.find(cutFaceGlobal[0]);
426                         bool otherSideFound = false;
428                         if (mpfAddrIter != masterPointFaceAddr.end())
429                         {
430                             bool miss = false;
432                             // Create the label count using point zero
433                             const labelList& masterFacesOfPZero = mpfAddrIter();
435                             labelList hits(masterFacesOfPZero.size(), 1);
437                             for
438                             (
439                                 label pointI = 1;
440                                 pointI < cutFaceGlobal.size();
441                                 pointI++
442                             )
443                             {
444                                 Map<labelList>::const_iterator
445                                     mpfAddrPointIter =
446                                         masterPointFaceAddr.find
447                                         (
448                                             cutFaceGlobal[pointI]
449                                         );
451                                 if
452                                 (
453                                     mpfAddrPointIter
454                                  == masterPointFaceAddr.end()
455                                 )
456                                 {
457                                     // Point is off the master patch. Skip
458                                     miss = true;
459                                     break;
460                                 }
462                                 const labelList& curMasterFaces =
463                                     mpfAddrPointIter();
465                                 // For every current face, try to find it in the
466                                 // zero-list
467                                 forAll (curMasterFaces, i)
468                                 {
469                                     forAll (masterFacesOfPZero, j)
470                                     {
471                                         if
472                                         (
473                                             curMasterFaces[i]
474                                          == masterFacesOfPZero[j]
475                                         )
476                                         {
477                                             hits[j]++;
478                                             break;
479                                         }
480                                     }
481                                 }
482                             }
484                             // If all point are found attempt matching
485                             if (!miss)
486                             {
487                                 forAll (hits, pointI)
488                                 {
489                                     if (hits[pointI] == cutFaceGlobal.size())
490                                     {
491                                         // Found other side.
492                                         otherSideFound = true;
494                                         cfMaster.append
495                                         (
496                                             masterFacesOfPZero[pointI]
497                                         );
499                                         cfSlave.append(faceI);
501                                         // Reverse the face such that it
502                                         // points out of the master patch
503                                         cf[cf.size() - 1] =
504                                             cf[cf.size() - 1].reverseFace();
506                                         if (debug)
507                                         {
508                                             Pout<< " other side: "
509                                                 << masterFacesOfPZero[pointI]
510                                                 << endl;
511                                         }
512                                     } // end of hits
513                                 } // end of for all hits
515                             } // end of not miss
517                             if (!otherSideFound || miss)
518                             {
519                                 if (debug)
520                                 {
521                                     Pout << " solo slave A" << endl;
522                                 }
524                                 cfMaster.append(-1);
525                                 cfSlave.append(faceI);
526                             }
527                         }
528                         else
529                         {
530                             // First point not in master patch
531                             if (debug)
532                             {
533                                 Pout << " solo slave B" << endl;
534                             }
536                             cfMaster.append(-1);
537                             cfSlave.append(faceI);
538                         }
539                     }
540                     else
541                     {
542                         if (debug)
543                         {
544                             Pout << " master side" << endl;
545                         }
547                         cfMaster.append(faceI - slavePatch_.size());
548                         cfSlave.append(-1);
549                     }
550                 }
551                 else
552                 {
553                     // Face not completed.  Prepare for the next point search
554                     prevPointLabel = curPointLabel;
555                     curPointLabel = bestAtanPoint;
556                     curPoint = lp[curPointLabel];
557                     cutFaceGlobalPoints.append(mp[curPointLabel]);
558                     cutFaceLocalPoints.append(curPointLabel);
560                     if (debug || cutFaceGlobalPoints.size() > faceSizeDebug)
561                     {
562                         faceSizeDebug *= 2;
564                         // Check for duplicate points in the face
565                         forAll (cutFaceGlobalPoints, checkI)
566                         {
567                             for
568                             (
569                                 label checkJ = checkI + 1;
570                                 checkJ < cutFaceGlobalPoints.size();
571                                 checkJ++
572                             )
573                             {
574                                 if
575                                 (
576                                     cutFaceGlobalPoints[checkI]
577                                  == cutFaceGlobalPoints[checkJ]
578                                 )
579                                 {
580                                     // Shrink local points for debugging
581                                     cutFaceLocalPoints.shrink();
583                                     face origFace;
584                                     face origFaceLocal;
585                                     if (faceI < slavePatch_.size())
586                                     {
587                                         origFace = slavePatch_[faceI];
588                                         origFaceLocal =
589                                             slavePatch_.localFaces()[faceI];
590                                     }
591                                     else
592                                     {
593                                         origFace =
594                                             masterPatch_
595                                             [faceI - slavePatch_.size()];
597                                         origFaceLocal =
598                                             masterPatch_.localFaces()
599                                             [faceI - slavePatch_.size()];
600                                     }
602                                     FatalErrorIn
603                                     (
604                                         "void enrichedPatch::"
605                                         "calcCutFaces() const"
606                                     )   << "Duplicate point found in cut face. "
607                                         << cutFaceGlobalPoints[checkI]
608                                         << " and " << checkI
609                                         << ".  Error in the face cutting "
610                                         << "algorithm for global face "
611                                         << origFace << " local face "
612                                         << origFaceLocal << nl
613                                         << "Slave size: " << slavePatch_.size()
614                                         << " Master size: "
615                                         << masterPatch_.size()
616                                         << " index: " << faceI << ".\n"
617                                         << "Face: " << curGlobalFace << nl
618                                         << "Cut face: " << cutFaceGlobalPoints
619                                         << " local: " << cutFaceLocalPoints
620                                         << nl << "Points: "
621                                         << face(cutFaceLocalPoints).points(lp)
622                                         << abort(FatalError);
623                                 }
624                             }
625                         }
626                     } // end of debug
627                 }
628             } while (!completedCutFace);
629         } // end of used edges
631         if (debug)
632         {
633             Pout << " Finished face " << faceI << endl;
634         }
636     } // end of local faces
638     // Re-pack the list into compact storage
639     cutFacesPtr_ = new faceList();
640     cutFacesPtr_->transfer(cf.shrink());
642     cutFaceMasterPtr_ = new labelList();
643     cutFaceMasterPtr_->transfer(cfMaster.shrink());
645     cutFaceSlavePtr_ = new labelList();
646     cutFaceSlavePtr_->transfer(cfSlave.shrink());
650 void Foam::enrichedPatch::calcLocalCutFaces() const
652     if (localCutFacesPtr_)
653     {
654         FatalErrorIn("void enrichedPatch::calcLocalCutFaces() const")
655             << "Local faces already calculated."
656             << abort(FatalError);
657     }
659     // Invert mesh points and renumber faces using it
660     const labelList& mp = meshPoints();
662     Map<label> mpLookup(2*mp.size());
664     forAll (mp, mpI)
665     {
666         mpLookup.insert(mp[mpI], mpI);
667     }
669     const faceList& cFaces = cutFaces();
671     localCutFacesPtr_ = new faceList(cFaces.size());
672     faceList& lcf = *localCutFacesPtr_;
674     forAll (cFaces, faceI)
675     {
676         const face& cf = cFaces[faceI];
678         face& curlcf = lcf[faceI];
680         curlcf.setSize(cf.size());
682         forAll (cf, pointI)
683         {
684             curlcf[pointI] = mpLookup.find(cf[pointI])();
685         }
686     }
690 void Foam::enrichedPatch::clearCutFaces()
692     deleteDemandDrivenData(cutFacesPtr_);
694     deleteDemandDrivenData(cutFaceMasterPtr_);
695     deleteDemandDrivenData(cutFaceSlavePtr_);
696     deleteDemandDrivenData(localCutFacesPtr_);
700 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
702 const Foam::faceList& Foam::enrichedPatch::cutFaces() const
704     if (!cutFacesPtr_)
705     {
706         calcCutFaces();
707     }
709     return *cutFacesPtr_;
713 const Foam::labelList& Foam::enrichedPatch::cutFaceMaster() const
715     if (!cutFaceMasterPtr_)
716     {
717         calcCutFaces();
718     }
720     return *cutFaceMasterPtr_;
724 const Foam::labelList& Foam::enrichedPatch::cutFaceSlave() const
726     if (!cutFaceSlavePtr_)
727     {
728         calcCutFaces();
729     }
731     return *cutFaceSlavePtr_;
735 const Foam::faceList& Foam::enrichedPatch::localCutFaces() const
737     if (!localCutFacesPtr_)
738     {
739         calcLocalCutFaces();
740     }
742     return *localCutFacesPtr_;
746 // ************************************************************************* //