Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / OpenFOAM / fields / pointPatchFields / pointPatchField / pointConstraint / pointConstraintI.H
blobe650ccb60438b0743c55025d5435cfefcc7b210d
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2010 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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
28 inline Foam::pointConstraint::pointConstraint()
30     Tuple2<label, vector>(0, vector::zero)
34 inline Foam::pointConstraint::pointConstraint(Istream& is)
36     Tuple2<label, vector>(is)
40 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
42 void Foam::pointConstraint::applyConstraint(const vector& cd)
44     if (first() == 0)
45     {
46         first() = 1;
47         second() = cd;
48     }
49     else if (first() == 1)
50     {
51         vector planeNormal = cd ^ second();
52         scalar magPlaneNormal = mag(planeNormal);
54         if (magPlaneNormal > 1e-3)
55         {
56             first() = 2;
57             second() = planeNormal/magPlaneNormal;
58         }
59     }
60     else if (first() == 2)
61     {
62         if (mag(cd & second()) > 1e-3)
63         {
64             first() = 3;
65             second() = vector::zero;
66         }
67     }
71 void Foam::pointConstraint::combine(const pointConstraint& pc)
73     if (first() == 0)
74     {
75         operator=(pc);
76     }
77     else if (first() == 1)
78     {
79         // Save single normal
80         vector n = second();
81         // Apply to supplied point constaint
82         operator=(pc);
83         applyConstraint(n);
84     }
85     else if (first() == 2)
86     {
87         if (pc.first() == 0)
88         {}
89         else if (pc.first() == 1)
90         {
91             applyConstraint(pc.second());
92         }
93         else if (pc.first() == 2)
94         {
95             // Both constrained to line. Same (+-)direction?
96             if (mag(second() & pc.second()) <= (1.0-1e-3))
97             {
98                 // Different directions
99                 first() = 3;
100                 second() = vector::zero;
101             }
102         }
103         else
104         {
105             first() = 3;
106             second() = vector::zero;
107         }
108     }
112 Foam::tensor Foam::pointConstraint::constraintTransformation() const
114     if (first() == 0)
115     {
116         return I;
117     }
118     else if (first() == 1)
119     {
120         return I - sqr(second());
121     }
122     else if (first() == 2)
123     {
124         return sqr(second());
125     }
126     else
127     {
128         return tensor::zero;
129     }
133 // ************************************************************************* //