Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
Mesh_def.hpp
1#ifndef Mesh_def_hpp
2#define Mesh_def_hpp
3
4#include "Mesh_decl.hpp"
5
14
15namespace FEDD {
16
17template <class SC, class LO, class GO, class NO>
18Mesh<SC,LO,GO,NO>::Mesh():
19numElementsGlob_(0),
20mapUnique_(),
21mapRepeated_(),
22pointsRep_(),
23pointsUni_(),
24bcFlagRep_(),
25bcFlagUni_(),
26surfaceElements_(),
27elementMap_(),
28comm_(),
29pointsRepRef_(),
30pointsUniRef_(),
31AABBTree_()
32{
33
34 surfaceElements_.reset(new Elements());
35
36 elementsC_.reset(new Elements());
37
38 FEType_ = "P1"; // We generally assume the mesh to be p1. In case of P1 or Q2 the FEType is allways adjusted
39}
40
41template <class SC, class LO, class GO, class NO>
42Mesh<SC,LO,GO,NO>::Mesh(CommConstPtrConst_Type& comm):
43numElementsGlob_(0),
44mapUnique_(),
45mapRepeated_(),
46pointsRep_(),
47pointsUni_(),
48bcFlagRep_(),
49bcFlagUni_(),
50surfaceElements_(),
51elementMap_(),
52edgeMap_(),
53comm_(comm),
54pointsRepRef_(),
55pointsUniRef_(),
56AABBTree_()
57{
58 AABBTree_.reset(new AABBTree_Type());
59 surfaceElements_.reset(new Elements());
60
61 elementsC_.reset(new Elements());
62
63 FEType_ = "P1";
64}
65
66template <class SC, class LO, class GO, class NO>
68
69 ElementsPtr_Type elements = this->getElementsC();
70// this->elementFlag_.reset( new vec_int_Type( elements->numberElements(), 0 ) );
71 if (type == "TPM_square") {
72 double xRef, yRef;
74 for (int i=0; i<elements->numberElements(); i++) {
75 xRef = ( this->pointsRep_->at( elements->getElement(i).getNode(0) )[0] + this->pointsRep_->at( elements->getElement(i).getNode(1) )[0] + this->pointsRep_->at( elements->getElement(i).getNode(2) )[0] ) / 3.;
76 yRef = ( this->pointsRep_->at( elements->getElement(i).getNode(0) )[1] + this->pointsRep_->at( elements->getElement(i).getNode(1) )[1] + this->pointsRep_->at( elements->getElement(i).getNode(2) )[1] ) / 3.;
77 if ( xRef>=0.3 && xRef<=0.7) {
78 if ( yRef>= 0.6) {
79 elements->getElement(i).setFlag(1);
80// this->elementFlag_->at(i) = 1;
81 }
82 }
83 }
84 }
85 else if (type == "Excavation1"){
86 }
87 else{
88
89 }
90
91}
92
93template <class SC, class LO, class GO, class NO>
94vec_int_ptr_Type Mesh<SC,LO,GO,NO>::getElementsFlag() const{
95 TEUCHOS_TEST_FOR_EXCEPTION( true, std::runtime_error, "we are not using the correct flags here. use the flags of elementC_." );
96 vec_int_ptr_Type tmp;
97 return tmp;
98}
99
100template <class SC, class LO, class GO, class NO>
101typename Mesh<SC,LO,GO,NO>::MapConstPtr_Type Mesh<SC,LO,GO,NO>::getMapUnique() const{
102
103 return mapUnique_;
104}
106template <class SC, class LO, class GO, class NO>
107typename Mesh<SC,LO,GO,NO>::MapConstPtr_Type Mesh<SC,LO,GO,NO>::getMapRepeated() const{
108
109 return mapRepeated_;
110}
111
112template <class SC, class LO, class GO, class NO>
113typename Mesh<SC,LO,GO,NO>::MapConstPtr_Type Mesh<SC,LO,GO,NO>::getElementMap() const{
114 TEUCHOS_TEST_FOR_EXCEPTION( elementMap_.is_null(), std::runtime_error, "Element map of mesh does not exist." );
115 return elementMap_;
116}
118// edgeMap
119template <class SC, class LO, class GO, class NO>
120typename Mesh<SC,LO,GO,NO>::MapConstPtr_Type Mesh<SC,LO,GO,NO>::getEdgeMap(){
121 TEUCHOS_TEST_FOR_EXCEPTION( edgeMap_.is_null(), std::runtime_error, "Element map of mesh does not exist." );
122 return edgeMap_;
123}
124
125template <class SC, class LO, class GO, class NO>
126vec2D_dbl_ptr_Type Mesh<SC,LO,GO,NO>::getPointsRepeated() const{
127
128 return pointsRep_;
129}
131template <class SC, class LO, class GO, class NO>
132vec2D_dbl_ptr_Type Mesh<SC,LO,GO,NO>::getPointsUnique() const{
133
134 return pointsUni_;
135}
136
137template <class SC, class LO, class GO, class NO>
139
140 return bcFlagRep_;
141}
142
143template <class SC, class LO, class GO, class NO>
144vec_int_ptr_Type Mesh<SC,LO,GO,NO>::getBCFlagUnique() const{
145
146 return bcFlagUni_;
148
149template <class SC, class LO, class GO, class NO>
150typename Mesh<SC,LO,GO,NO>::ElementsPtr_Type Mesh<SC,LO,GO,NO>::getElementsC(){
151 return elementsC_;
152}
153
154template <class SC, class LO, class GO, class NO>
155typename Mesh<SC,LO,GO,NO>::ElementsPtr_Type Mesh<SC,LO,GO,NO>::getSurfaceElements(){
156 return surfaceElements_;
157}
158
159template <class SC, class LO, class GO, class NO>
161
162 return dim_;
163}
164
165template <class SC, class LO, class GO, class NO>
167
168 return numElementsGlob_;
169}
170
171template <class SC, class LO, class GO, class NO>
173 TEUCHOS_TEST_FOR_EXCEPTION( this->elementsC_.is_null(), std::runtime_error ,"Elements do not exist." );
174 return this->elementsC_->numberElements();
175}
176
177template <class SC, class LO, class GO, class NO>
179 if (!type.compare("Unique"))
180 return pointsUni_->size();
181 else if(!type.compare("Repeated"))
182 return pointsRep_->size();
183
184 TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,"Select valid map type: unique or repeated.");
185 return 0;
186}
188template <class SC, class LO, class GO, class NO>
191 switch (dim_) {
192 case 2:
193 if ( !FEType_.compare("P1") )
194 return 3;
195 else if ( !FEType_.compare("P1-disc") || !FEType_.compare("P1-disc-global") ){
196 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, "P1-disc only available in 3D.");
197 }
198 else if( !FEType_.compare("P2") )
199 return 6;
200 else if( !FEType_.compare("P2-CR") ){
201 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, "P2-CR only available in 3D.");
202 }
203 else if( !FEType_.compare("Q2-20") ){
204 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, "Q2-20 only available in 3D.");
205 }
206 else if( !FEType_.compare("Q2") ){
207 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, "Q2 only available in 3D.");
208 }
209 break;
210 case 3:
211 if ( !FEType_.compare("P1") )
212 return 4;
213 if ( !FEType_.compare("P1-disc") || !FEType_.compare("P1-disc-global") )
214 return 4;
215 else if( !FEType_.compare("P2") )
216 return 10;
217 else if( !FEType_.compare("P2-CR") )
218 return 15;
219 else if( !FEType_.compare("Q2-20") )
220 return 20; // Q2 Serendipity (missing interior volume node and missing interior face nodes gives 20 nodes)
221 else if( !FEType_.compare("Q2") )
222 return 27;
223 break;
224 default:
225 return -1;
226 break;
227 }
228 return -1;
229}
230
231template <class SC, class LO, class GO, class NO>
233{
234 // Bemerkung: Repeated und Unique sind unterschiedlich lang!!! => zwei Schleifen
235
236 // Setze zunaechst alles auf Null, andernfalls kann man nicht drauf zugreifen
237 // vec2D_dbl_ptr_Type zeroRep(new vec2D_dbl_Type(pointsRep_->size(),vec_dbl_Type(pointsRep_->at(0).size(),0.0)));
238 // vec2D_dbl_ptr_Type zeroUni(new vec2D_dbl_Type(pointsUni_->size(),vec_dbl_Type(pointsUni_->at(0).size(),0.0)));
239
240
241 pointsRepRef_.reset( new vec2D_dbl_Type() );
242 pointsUniRef_.reset( new vec2D_dbl_Type() );
243 // zeroRep und zeroUni leben nur hier drinnen, weswegen wir Pointer gleich Pointer setzen koennen.
244 // TODO: *PointsRepRef_ = *zeroRep funktioniert nicht.
245 *pointsRepRef_ = *pointsRep_;
246 *pointsUniRef_ = *pointsUni_;
247
248 // // Repeated
249 // for(int i = 0; i < pointsRep_->size(); i++)
250 // {
251 // for(int j = 0; j < pointsRep_->at(0).size(); j++)
252 // {
253 // pointsRepRef_->at(i).at(j) = PointsRep_->at(i).at(j);
254 // }
255 // }
256 //
257 // // Unique
258 // for(int i = 0; i < PointsUni_->size(); i++)
259 // {
260 // for(int j = 0; j < PointsUni_->at(0).size(); j++)
261 // {
262 // PointsUniRef_->at(i).at(j) = PointsUni_->at(i).at(j);
263 // }
264 // }
265}
266
267template <class SC, class LO, class GO, class NO>
268void Mesh<SC,LO,GO,NO>::moveMesh( MultiVectorPtr_Type displacementUnique, MultiVectorPtr_Type displacementRepeated )
269{
270 // Bemerkung: Repeated und Unique sind unterschiedlich lang!!! => zwei Schleifen
271 TEUCHOS_TEST_FOR_EXCEPTION (displacementRepeated.is_null(), std::runtime_error," displacementRepeated in moveMesh is null.")
272 TEUCHOS_TEST_FOR_EXCEPTION (displacementUnique.is_null(), std::runtime_error," displacementRepeated in moveMesh is null.")
273 // Repeated
274 Teuchos::ArrayRCP<const SC> values = displacementRepeated->getData(0); //only 1 MV
275 for(int i = 0; i < pointsRepRef_->size(); i++)
276 {
277 for(int j = 0; j < pointsRepRef_->at(0).size(); j++)
278 {
279 // Sortierung von DisplacementRepeated ist x-y-x-y-x-y-x-y bzw. x-y-z-x-y-z-x-y-z
280 // Achtung: DisplacementRepeated ist ein Pointer der mit (*) dereferenziert werden muss.
281 // Operator[] kann nicht auf einen Pointer angewendet werden!!!
282 // Es sei denn es ist ein Array.
283 pointsRep_->at(i).at(j) = pointsRepRef_->at(i).at(j) + values[dim_*i+j];
284 }
285 }
286
287 // Unique
288 values = displacementUnique->getData(0); //only 1 MV
289 for(int i = 0; i < pointsUniRef_->size(); i++)
290 {
291 for(int j = 0; j < pointsUniRef_->at(0).size(); j++)
292 {
293 // Sortierung von DisplacementRepeated ist x-y-x-y-x-y-x-y bzw. x-y-z-x-y-z-x-y-z
294 // Erklaerung: DisplacementUnique ist ein Vector-Pointer, wo in jedem Eintrag ein MultiVector-Pointer drin steht (std::vector<MultiVector_ptr_Type>)
295 // Greife mit ->at auf den Eintrag des Vektors zu (hier nur ein Eintrag vorhanden), dereferenziere den damit erhaltenen MultiVector-Pointer (als Referenz) um einen
296 // MultiVector zu erhalten.
297 // Greife dann mit [] auf das entsprechende Array (double *&) im MultiVector zu (hier gibt es nur einen)
298 // und anschliessend mit [] auf den Wert des Arrays.
299 // Beachte falls x ein Array ist (also z.B. double *), dann ist x[i] := *(x+i)!!!
300 // Liefert also direkt den Wert und keinen Pointer auf einen double.
301 // Achtung: MultiVector[] liefert double* wohingegen MultiVector() Epetra_Vector* zurueck liefert
302 pointsUni_->at(i).at(j) = pointsUniRef_->at(i).at(j) + values[dim_*i+j];
303 }
304 }
305}
306
307template <class SC, class LO, class GO, class NO>
308void Mesh<SC,LO,GO,NO>::create_AABBTree(){
309 if (AABBTree_.is_null()){
310 AABBTree_.reset(new AABBTree_Type());
311 }
312 AABBTree_->createTreeFromElements(
313 getElementsC(),
314 getPointsRepeated()
315 );
316}
317
318template <class SC, class LO, class GO, class NO>
319vec_int_ptr_Type Mesh<SC,LO,GO,NO>::findElemsForPoints(
320 vec2D_dbl_ptr_Type queryPoints
321){
322 int numPoints = queryPoints->size();
323
324 // Return vector. -1 means that point is in no elem, otherwise entry is the
325 // elem the point is in
326 vec_int_ptr_Type pointToElem(
327 new vec_int_Type(
328 numPoints,
329 -1
330 )
331 );
332
333 // Create tree if it is empty
334 if (AABBTree_->isEmpty()){
335 AABBTree_->createTreeFromElements(
336 getElementsC(),
337 getPointsRepeated(),
338 false
339 );
340 }
341
342 // Query the AABBTree
343 std::map<int, std::list<int> > treeToItem;
344 std::map<int, std::list<int> > itemToTree;
345 tie(treeToItem, itemToTree) = AABBTree_->scanTree(queryPoints, false);
346
347 // FIXME: put this in a function of AABBTree?
348 // unnest the returned answer for each query_point
349 int point = -1;
350 bool found = false;
351 std::list<int> rectangles;
352 std::list<int> elements;
353 for (auto keyValue: itemToTree){
354 // FIXME: put this in a function of AABBTree?
355 // rectangles is a list<int> of all rectangles point is in
356 // find the element(s) that is/are in all of said rectangles.
357 // If there is only one element that is the element the point is in,
358 // if not we have to query all remaining elements
359 point = keyValue.first;
360 rectangles = keyValue.second;
361
362
363 // query all remaining elements
364 for(auto rectangle: rectangles){
365 elements = AABBTree_->getElements(rectangle);
366 for(auto element: elements){
367 found = isPointInElem(queryPoints->at(point), element);
368 if( found ){
369 pointToElem->at(point) = element;
370 break;
371 }
372 }
373 if( found ){
374 // we already found the element, no need to check additional rectangles
375 break;
376 }
377 }
378
379 }
380 return pointToElem;
381}
382
383template <class SC, class LO, class GO, class NO>
384vec_dbl_Type Mesh<SC,LO,GO,NO>::getBaryCoords(vec_dbl_Type point, int element){
385 vec_int_Type localNodes = elementsC_->getElement(element).getVectorNodeList();
386 vec2D_dbl_Type coords(
387 localNodes.size(),
388 vec_dbl_Type(2, 0.0) // FIXME: this depends on the dimension
389 );
390 int node = 0;
391 for (int localNode=0; localNode<localNodes.size(); localNode++){
392 node = localNodes.at(localNode); // get global node_id
393 coords.at(localNode) = pointsRep_->at(node); //get global coordinates
394 }
395
396 double px, py, x1, x2, x3, y1, y2, y3;
397 px = point.at(0); py = point.at(1);
398 x1 = coords.at(0).at(0); y1 = coords.at(0).at(1);
399 x2 = coords.at(1).at(0); y2 = coords.at(1).at(1);
400 x3 = coords.at(2).at(0); y3 = coords.at(2).at(1);
401
402 // baryzentric coordinates
403 double det_T = (y2 - y3) * (x1 - x3) + (x3 - x2) * (y1 - y3);
404
405 vec_dbl_Type baryCoords(
406 3,
407 0.0
408 );
409 baryCoords[0] =(y2 - y3) * (px - x3) + (x3 - x2) * (py - y3);
410 baryCoords[0] = baryCoords[0] / det_T;
411
412 baryCoords[1] = (y3 -y1) * (px - x3) + (x1 - x3) * (py - y3);
413 baryCoords[1] = baryCoords[1] / det_T;
414
415 baryCoords[2] = 1 - baryCoords[1] - baryCoords[0];
416 return baryCoords;
417}
418
419
420template <class SC, class LO, class GO, class NO>
421bool Mesh<SC,LO,GO,NO>:: isPointInElem(vec_dbl_Type point, int element){
422 // FIXME: This is only valid for a triangle
423 vec_dbl_Type baryCoords;
424 baryCoords = getBaryCoords(point, element);
425
426 if(baryCoords[0] >= 0 && baryCoords[1] >= 0 && baryCoords[2] >= 0){
427 return true;
428 }
429 return false;
430}
431
432template <class SC, class LO, class GO, class NO>
434 this->elementsVec_ = Teuchos::rcp( new vec2D_int_Type( this->elementsC_->numberElements() ) );
435 for (int i=0; i<this->elementsVec_->size(); i++)
436 this->elementsVec_->at(i) = this->elementsC_->getElement(i).getVectorNodeList();
437
438 return this->elementsVec_;
439}
440
441template <class SC, class LO, class GO, class NO>
443
444 int outwardNormals = 0;
445 int inwardNormals = 0;
446 for (UN T=0; T<elementsC_->numberElements(); T++) {
447 FiniteElement fe = elementsC_->getElement( T );
448 ElementsPtr_Type subEl = fe.getSubElements(); // might be null
449 for (int surface=0; surface<fe.numSubElements(); surface++) {
450 FiniteElement feSub = subEl->getElement( surface );
451 vec_int_Type nodeListElement = fe.getVectorNodeList();
452 if(subEl->getDimension() == dim_-1 ){
453 vec_int_Type nodeList = feSub.getVectorNodeListNonConst();
454 int numNodes_T = nodeList.size();
455
456 vec_dbl_Type v_E(dim_,1.);
457 double norm_v_E=1.;
458 LO id0 = nodeList[0];
459
460 Helper::computeSurfaceNormal(dim_, pointsRep_,nodeList,v_E,norm_v_E);
461
462 std::sort(nodeList.begin(), nodeList.end());
463 std::sort(nodeListElement.begin(), nodeListElement.end());
464
465 std::vector<int> v_symDifference;
466
467 std::set_symmetric_difference(
468 nodeList.begin(), nodeList.end(),
469 nodeListElement.begin(), nodeListElement.end(),
470 std::back_inserter(v_symDifference));
471
472 LO id1 = v_symDifference[0]; // This is a node that is not part of the surface, i.e. 4th element node in 3D
473
474 vec_dbl_Type p0(dim_,0.);
475 for(int i=0; i< dim_; i++)
476 p0[i] = pointsRep_->at(id1)[i] - pointsRep_->at(id0)[i];
477
478 double sum = 0.;
479 for(int i=0; i< dim_; i++)
480 sum += p0[i] * v_E[i];
481
482 if(sum<=0){
483 outwardNormals++;
484 }
485 if(sum>0){
486 inwardNormals++;
487 }
488 if(sum>0) // if the sum is greater than 0, the normal is in inward direction and thus is flipped.
489 flipSurface(subEl,surface);
490
491
492 }
493 }
494 }
495 Teuchos::reduceAll<int, int> (*this->getComm(), Teuchos::REDUCE_SUM, inwardNormals, Teuchos::outArg (inwardNormals));
496 Teuchos::reduceAll<int, int> (*this->getComm(), Teuchos::REDUCE_SUM, outwardNormals, Teuchos::outArg (outwardNormals));
497
498 if(this->getComm()->getRank() == 0){
499 std::cout << " ############################################ " << std::endl;
500 std::cout << " Mesh Orientation Statistic " << std::endl;
501 std::cout << " Number of outward normals " << outwardNormals << std::endl;
502 std::cout << " Number of inward normals " << inwardNormals << std::endl;
503 std::cout << " ############################################ " << std::endl;
504 }
505
506}
507template <class SC, class LO, class GO, class NO>
509
510 int posDet=0;
511 int negDet=0;
512
513 SC detB;
514 SmallMatrix<SC> B(dim_);
515 SmallMatrix<SC> Binv(dim_);
516
517 for (UN T=0; T<elementsC_->numberElements(); T++) {
518 Helper::buildTransformation(elementsC_->getElement(T).getVectorNodeList(), pointsRep_, B);
519 detB = B.computeInverse(Binv);
520
521 if(detB<0){
522 negDet++;
523 }
524 if(detB>0){
525 posDet++;
526 }
527 if(detB<0) // If the determinant is smaller than zero we flip the element
528 flipElement(elementsC_,T);
529
530 }
531 std::cout << " Finished " << std::endl;
532 Teuchos::reduceAll<int, int> (*this->getComm(), Teuchos::REDUCE_SUM, negDet, Teuchos::outArg (negDet));
533 Teuchos::reduceAll<int, int> (*this->getComm(), Teuchos::REDUCE_SUM, posDet, Teuchos::outArg (posDet));
534
535 if(this->getComm()->getRank() == 0){
536 std::cout << " ############################################ " << std::endl;
537 std::cout << " Mesh Orientation Statistic " << std::endl;
538 std::cout << " Number of positive dets " << posDet << std::endl;
539 std::cout << " Number of negative dets " << negDet << std::endl;
540 std::cout << " ############################################ " << std::endl;
541 }
542
543}
544
545// We allways want a outward normal direction
546// Assumptions: We are flipping a surface which is a subelement. Subelement are generally element that are on the boundary layers of the domain. Thus, they are unique. (There are no two identical triangles in two elements as subelements)
547// Question: Easiest way to flip the surface without redoing whole dim-element (surface being dim-1-element)
548template <class SC, class LO, class GO, class NO>
549void Mesh<SC,LO,GO,NO>::flipSurface(ElementsPtr_Type subEl, int surfaceNumber){
550
551 vec_LO_Type surfaceElements_vec = subEl->getElement(surfaceNumber).getVectorNodeList();
552
553 if(dim_ == 2){
554 if(FEType_ == "P1"){
555 LO id1,id2;
556 id1= surfaceElements_vec[0];
557 id2= surfaceElements_vec[1];
558
559 surfaceElements_vec[0] = id2;
560 surfaceElements_vec[1] = id1;
561 }
562 else if(FEType_ == "P2"){
563 LO id1,id2,id3;
564 id1= surfaceElements_vec[0];
565 id2= surfaceElements_vec[1];
566 id3= surfaceElements_vec[2];
567
568 surfaceElements_vec[0] = id2;
569 surfaceElements_vec[1] = id1;
570 surfaceElements_vec[2] = id3;
571 }
572 }
573 else if(dim_ == 3){
574
575 if(FEType_ == "P1"){
576 LO id1,id2,id3,id4,id5,id6;
577 id1= surfaceElements_vec[0];
578 id2= surfaceElements_vec[1];
579 id3= surfaceElements_vec[2];
580
581 surfaceElements_vec[0] = id1;
582 surfaceElements_vec[1] = id3;
583 surfaceElements_vec[2] = id2;
584 }
585 else if(FEType_ == "P2"){
586 LO id1,id2,id3,id4,id5,id6;
587 id1= surfaceElements_vec[0];
588 id2= surfaceElements_vec[1];
589 id3= surfaceElements_vec[2];
590 id4= surfaceElements_vec[3];
591 id5= surfaceElements_vec[4];
592 id6= surfaceElements_vec[5];
593
594 surfaceElements_vec[0] = id1;
595 surfaceElements_vec[1] = id3;
596 surfaceElements_vec[2] = id2;
597 surfaceElements_vec[3] = id6;
598 surfaceElements_vec[4] = id5;
599 surfaceElements_vec[5] = id4;
600 }
601 else
602 TEUCHOS_TEST_FOR_EXCEPTION( true, std::runtime_error, "We can only flip normals for P1 or P2 elements. Invalid " << FEType_ << " " );
603
604 }
605 FiniteElement feFlipped(surfaceElements_vec,subEl->getElement(surfaceNumber).getFlag());
606 subEl->switchElement(surfaceNumber,feFlipped); // We can switch the current element with the newly defined element which has just a different node ordering.
607
608
609}
610
611// We allways want a positive determinant
612template <class SC, class LO, class GO, class NO>
613void Mesh<SC,LO,GO,NO>::flipElement(ElementsPtr_Type elements, int elementNumber){
614
615 vec_LO_Type surfaceElements_vec = elements->getElement(elementNumber).getVectorNodeList();
616 if(dim_ == 2){
617 if(FEType_ == "P1"){
618 LO id1,id2,id3;
619 id1= surfaceElements_vec[0];
620 id2= surfaceElements_vec[1];
621 id3= surfaceElements_vec[2];
622
623 surfaceElements_vec[0] = id1;
624 surfaceElements_vec[1] = id3;
625 surfaceElements_vec[2] = id2;
626
627 }
628 else if(FEType_ == "P2"){
629 LO id1,id2,id3,id4,id5,id6;
630 id1= surfaceElements_vec[0];
631 id2= surfaceElements_vec[1];
632 id3= surfaceElements_vec[2];
633 id4= surfaceElements_vec[3];
634 id5= surfaceElements_vec[4];
635 id6= surfaceElements_vec[5];
636
637 surfaceElements_vec[0] = id1;
638 surfaceElements_vec[1] = id3;
639 surfaceElements_vec[2] = id2;
640 surfaceElements_vec[3] = id6;
641 surfaceElements_vec[4] = id5;
642 surfaceElements_vec[5] = id4;
643 }
644 }
645 else if(dim_ == 3){
646
647 if(FEType_ == "P1"){
648 LO id1,id2,id3,id4;
649 id1= surfaceElements_vec[0];
650 id2= surfaceElements_vec[1];
651 id3= surfaceElements_vec[2];
652 id4= surfaceElements_vec[3];
653
654 surfaceElements_vec[0] = id1;
655 surfaceElements_vec[1] = id3;
656 surfaceElements_vec[2] = id2;
657 surfaceElements_vec[3] = id4;
658
659 }
660 else if(FEType_ == "P2"){
661 LO id1,id2,id3,id4,id5,id6, id7,id8,id9,id0;
662 id0= surfaceElements_vec[0];
663 id1= surfaceElements_vec[1];
664 id2= surfaceElements_vec[2];
665 id3= surfaceElements_vec[3];
666 id4= surfaceElements_vec[4];
667 id5= surfaceElements_vec[5];
668 id6= surfaceElements_vec[6];
669 id7= surfaceElements_vec[7];
670 id8= surfaceElements_vec[8];
671 id9= surfaceElements_vec[9];
672
673 surfaceElements_vec[0] = id0;
674 surfaceElements_vec[1] = id2;
675 surfaceElements_vec[2] = id1;
676 surfaceElements_vec[3] = id3;
677
678 surfaceElements_vec[4] = id6;
679 surfaceElements_vec[5] = id5;
680 surfaceElements_vec[6] = id4;
681 surfaceElements_vec[7] = id7;
682 surfaceElements_vec[8] = id9;
683 surfaceElements_vec[9] = id8;
684 }
685 else
686 TEUCHOS_TEST_FOR_EXCEPTION( true, std::runtime_error, "We can only flip normals for P1 or P2 elements. Invalid " << FEType_ << " " );
687
688 }
689 FiniteElement feFlipped(surfaceElements_vec,elements->getElement(elementNumber).getFlag());
690 ElementsPtr_Type subElements = elements->getElement(elementNumber).getSubElements();
691 for(int T = 0; T<elements->getElement(elementNumber).numSubElements(); T++){
692 if(!feFlipped.subElementsInitialized())
693 feFlipped.initializeSubElements( this->FEType_, this->dim_ -1) ;
694 feFlipped.addSubElement(subElements->getElement(T));
695 }
696 elements->switchElement(elementNumber,feFlipped); // We can switch the current element with the newly defined element which has just a different node ordering.
697
698}
699
700}
701#endif
Definition FiniteElement.hpp:17
static void buildTransformation(const vec_int_Type &element, vec2D_dbl_ptr_Type pointsRep, SmallMatrix< SC > &B, std::string FEType="P")
Build transformation of element to reference element depending on FEType.
Definition Helper.cpp:138
static void computeSurfaceNormal(int dim, vec2D_dbl_ptr_Type pointsRep, vec_int_Type nodeList, vec_dbl_Type &v_E, double &norm_v_E)
Compute surface normal of corresponding surface.
Definition Helper.cpp:103
ElementsPtr_Type getElementsC()
Returns element list as c-object.
Definition Mesh_def.hpp:150
vec2D_dbl_ptr_Type getPointsRepeated() const
getter for list of repeated points with x,y,z coordinates in each row
Definition Mesh_def.hpp:126
CommConstPtrConst_Type getComm()
Communicator object.
Definition Mesh_decl.hpp:138
MapConstPtr_Type getEdgeMap()
Getter for edge map.
Definition Mesh_def.hpp:120
vec2D_dbl_ptr_Type getPointsUnique() const
getter for list of unique points with x,y,z coordinates in each row
Definition Mesh_def.hpp:132
vec_int_ptr_Type getBCFlagRepeated() const
Getter for flags corresponting to repeated points.
Definition Mesh_def.hpp:138
void setElementFlags(std::string type="")
Something for TPM. Can be deprecated soon.
Definition Mesh_def.hpp:67
void correctNormalDirections()
Correcting the normal direction of all surface normals set as subelements of volume elements to be ou...
Definition Mesh_def.hpp:442
GO getNumElementsGlobal()
Global number of elements.
Definition Mesh_def.hpp:166
vec_int_ptr_Type getElementsFlag() const
Getter for element flags.
Definition Mesh_def.hpp:94
MapConstPtr_Type getMapUnique() const
Getter for unique node map.
Definition Mesh_def.hpp:101
MapConstPtr_Type getElementMap() const
Getter for element map.
Definition Mesh_def.hpp:113
void correctElementOrientation()
Correct the element orientation of all elements to have positive volume / det when doint transformati...
Definition Mesh_def.hpp:508
vec2D_int_ptr_Type getElements()
Returns elements as a vector type contrary to the C-object list.
Definition Mesh_def.hpp:433
ElementsPtr_Type getSurfaceElements()
Getter for surface elements. Probably set in mesh partitioner. They are generally the dim-1 surface e...
Definition Mesh_def.hpp:155
MapConstPtr_Type getMapRepeated() const
Getter for repeated node mal.
Definition Mesh_def.hpp:107
LO getNumElements()
Local number of elements.
Definition Mesh_def.hpp:172
LO getNumPoints(std::string type="Unique")
Get local (LO) number of points either in unique or repeated version.
Definition Mesh_def.hpp:178
int getDimension()
Definition Mesh_def.hpp:160
void moveMesh(MultiVectorPtr_Type displacementUnique, MultiVectorPtr_Type displacementRepeated)
Moving mesh according to displacement based on reference configuration.
Definition Mesh_def.hpp:268
int getOrderElement()
Definition Mesh_def.hpp:189
void setReferenceConfiguration()
Setting current coordinates as reference configuration. Should only be called once :D.
Definition Mesh_def.hpp:232
vec_int_ptr_Type getBCFlagUnique() const
Getter for flags corresponting to unique points.
Definition Mesh_def.hpp:144
This class represents a templated small Matrix of type T. Primarily created for 2x2 and 3x3 matrices....
Definition SmallMatrix.hpp:22
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement.cpp:5