Salome HOME
Class GenericMatrix no longer uses DoubleTab
[tools/solverlab.git] / CDMATH / mesh / src / Node.cxx
1 /*
2  * Node.cxx
3  *
4  *  Created on: 23 janv. 2012
5  *      Authors: CDMAT
6  */
7
8 #include "Node.hxx"
9 #include "CdmathException.hxx"
10 #include <algorithm> 
11
12 //----------------------------------------------------------------------
13 Node::Node( void )
14 //----------------------------------------------------------------------
15 {
16         _numberOfCells = 0 ;
17         _numberOfFaces = 0 ;
18         _numberOfEdges = 0 ;
19         _groupNames=std::vector<std::string>(0);
20         _region=-1;
21     _isBorder=false;
22 }
23
24 //----------------------------------------------------------------------
25 Node::~Node( void )
26 //----------------------------------------------------------------------
27 {
28 }
29
30 //----------------------------------------------------------------------
31 Node::Node( const Node& node )
32 //----------------------------------------------------------------------
33 {
34         _point = node.getPoint() ;
35         _cellsId = node.getCellsId() ;
36         _facesId = node.getFacesId() ;
37         _numberOfCells = node.getNumberOfCells() ;
38         _numberOfFaces = node.getNumberOfFaces() ;
39         _numberOfEdges = node.getNumberOfEdges() ;
40         _groupNames=node.getGroupNames();
41         _region=node.getRegion();
42     _isBorder=node.isBorder();
43 }
44
45 //----------------------------------------------------------------------
46 Node::Node( const int numberOfCells, const int numberOfFaces, const int numberOfEdges, const Point p)
47 //----------------------------------------------------------------------
48 {
49
50         _point = p ;
51         _numberOfCells = numberOfCells ;
52         _numberOfFaces = numberOfFaces ;
53         _numberOfEdges = numberOfEdges ;
54         _cellsId = std::vector< int >(_numberOfCells,0);
55         _facesId = std::vector< int >(_numberOfFaces,0);
56         _neighbourNodesId = std::vector< int >(_numberOfEdges,0);
57         _groupNames=std::vector<std::string>(0);
58         _region=-1;
59     _isBorder=false;
60 }
61
62 //----------------------------------------------------------------------
63 std::vector< int >
64 Node::getCellsId( void ) const 
65 //----------------------------------------------------------------------
66 {
67         return _cellsId ;
68 }
69
70 //----------------------------------------------------------------------
71 std::vector< int >
72 Node::getFacesId( void ) const 
73 //----------------------------------------------------------------------
74 {
75         return _facesId ;
76 }
77
78 //----------------------------------------------------------------------
79 std::vector< int >
80 Node::getNeighbourNodesId( void ) const 
81 //----------------------------------------------------------------------
82 {
83         return _neighbourNodesId ;
84 }
85
86 //----------------------------------------------------------------------
87 int
88 Node::getNumberOfCells( void ) const 
89 //----------------------------------------------------------------------
90 {
91         return _numberOfCells ;
92 }
93
94 //----------------------------------------------------------------------
95 int
96 Node::getNumberOfFaces( void ) const 
97 //----------------------------------------------------------------------
98 {
99         return _numberOfFaces ;
100 }
101
102 //----------------------------------------------------------------------
103 int
104 Node::getNumberOfEdges( void ) const 
105 //----------------------------------------------------------------------
106 {
107         return _numberOfEdges ;
108 }
109
110 //----------------------------------------------------------------------
111 Point
112 Node::getPoint( void ) const
113 //----------------------------------------------------------------------
114 {
115         return _point ;
116 }
117
118 //----------------------------------------------------------------------
119 double
120 Node::x( void ) const 
121 //----------------------------------------------------------------------
122 {
123         return _point.x() ;
124 }
125
126 //----------------------------------------------------------------------
127 double
128 Node::y( void ) const 
129 //----------------------------------------------------------------------
130 {
131         return _point.y() ;
132 }
133
134 //----------------------------------------------------------------------
135 double
136 Node::z( void ) const 
137 //----------------------------------------------------------------------
138 {
139         return _point.z() ;
140 }
141
142 std::vector<std::string>
143 Node::getGroupNames(void) const
144 {
145         return _groupNames;
146 }
147
148 std::string
149 Node::getGroupName(int igroup) const
150 {
151     if (igroup<_groupNames.size())
152         return _groupNames[igroup];
153     else
154     {
155         std::cout<<"Error Node::getGroupName(int igroup), group number "<<igroup+1<<" requested"<<std::endl;
156         std::cout<<"Node belongs to "<< _groupNames.size() <<" groups"<<std::endl;
157         throw CdmathException("Node has no group with number igroup");
158     }
159 }
160
161 void
162 Node::setGroupName(const std::string groupName)
163 {
164         _groupNames.insert(_groupNames.begin(),groupName);
165         _region=0;
166 }
167
168 bool
169 Node::isBorder(void) const
170 {
171         if (_region==0 | _isBorder)
172                 return true;
173         else
174                 return false;
175 }
176
177 int
178 Node::getRegion(void) const
179 {
180         return _region;
181 }
182
183 //----------------------------------------------------------------------
184 void
185 Node::addFaceId (const int numFace, const int faceId, bool isBorder  )
186 //----------------------------------------------------------------------
187 {
188         _facesId[numFace] = faceId ;
189     if(isBorder)
190     {
191         if(std::find(_groupNames.begin(), _groupNames.end(), "Boundary") == _groupNames.end())//No group named Boundary
192             _groupNames.insert(_groupNames.begin(),"Boundary");
193         _isBorder=true;
194     }
195 }
196
197 //----------------------------------------------------------------------
198 void
199 Node::addCellId (const int numCell, const int cellId )
200 //----------------------------------------------------------------------
201 {
202         _cellsId[numCell] = cellId ;
203 }
204
205 //----------------------------------------------------------------------
206 void
207 Node::addNeighbourNodeId (const int numNode, const int nodeId ) 
208 //----------------------------------------------------------------------
209 {
210         _neighbourNodesId[numNode] = nodeId ;
211 }
212
213 //----------------------------------------------------------------------
214 double
215 Node::distance( const Node& n ) const
216 //----------------------------------------------------------------------
217 {
218         double distance=_point.distance(n.getPoint());
219         return distance ;
220 }
221
222 //----------------------------------------------------------------------
223 int
224 //----------------------------------------------------------------------
225 Node::getFaceId(int localId) const
226 //----------------------------------------------------------------------
227 {
228     if(localId<_numberOfFaces)
229         return _facesId[localId];
230     else
231     {
232         std::cout<< "Local id requested : "<< localId<<", total number of faces= "<<_numberOfFaces<<std::endl;
233         throw CdmathException("Node::getFaceId : incorrect face local id");
234     }
235 }
236
237 //----------------------------------------------------------------------
238 int
239 //----------------------------------------------------------------------
240 Node::getNeighbourNodeId(int localId) const
241 //----------------------------------------------------------------------
242 {
243     if(localId<_numberOfEdges)
244         return _neighbourNodesId[localId];
245     else
246     {
247         std::cout<< "Local id requested : "<< localId<<", total number of edges= "<<_numberOfEdges<<std::endl;
248         throw CdmathException("Node::getNeighbourNodeId : incorrect node local id");
249     }
250 }
251
252 const Node& 
253 Node::operator= ( const Node& node )
254 //----------------------------------------------------------------------
255 {
256    _point = node.getPoint()  ;
257    _cellsId = node.getCellsId() ;
258    _facesId = node.getFacesId() ;
259    _neighbourNodesId = node.getNeighbourNodesId() ;
260    _numberOfCells = node.getNumberOfCells() ;
261    _numberOfFaces = node.getNumberOfFaces() ;
262    _numberOfEdges = node.getNumberOfEdges() ;
263    _groupNames = node.getGroupNames();  
264    _isBorder=node.isBorder();
265         return *this;
266 }