LR-splines  0.5
MeshRectangle.cpp
Go to the documentation of this file.
00001 
00002 #include "LRSpline/MeshRectangle.h"
00003 #include "LRSpline/Element.h"
00004 #include "LRSpline/Basisfunction.h"
00005 #include <algorithm>
00006 
00007 namespace LR {
00008 
00009 #define DOUBLE_TOL 1e-14
00010 
00011 bool MeshRectangle::addUniqueRect(std::vector<MeshRectangle*> &rects, MeshRectangle* newRect) {
00012         for(MeshRectangle *m : rects) {
00013                 if(m->equals(newRect) ) {
00014                         return false;
00015                 }
00016         }
00017         rects.push_back(newRect);
00018         return true;
00019 };
00020 
00021 MeshRectangle::MeshRectangle() {
00022         start_.resize(3);
00023         stop_.resize(3);
00024         start_[0]     = 0;
00025         start_[1]     = 0;
00026         start_[2]     = 0;
00027         stop_[0]      = 0;
00028         stop_[1]      = 0;
00029         stop_[2]      = 0;
00030         multiplicity_ = 0;
00031         constDir_     = -1;
00032 }
00033 
00034 MeshRectangle::MeshRectangle(double start_u,
00035                              double start_v,
00036                              double start_w,
00037                              double stop_u,
00038                              double stop_v,
00039                              double stop_w,
00040                              int multiplicity) {
00041         start_.resize(3);
00042         stop_.resize(3);
00043         start_[0]     = start_u;
00044         start_[1]     = start_v;
00045         start_[2]     = start_w;
00046         stop_[0]      = stop_u; 
00047         stop_[1]      = stop_v; 
00048         stop_[2]      = stop_w; 
00049         multiplicity_ =  multiplicity;
00050         constDir_     = -1;
00051 
00052         for(int i=0; i<3; i++)
00053                 if(start_[i] == stop_[i])
00054                         constDir_ = i;
00055         if(constDir_ == -1)
00056                 std::cerr << "Error creating MeshRectangle: Not parallel to the parametric axis\n";
00057 }
00058 
00059 MeshRectangle::~MeshRectangle() {
00060 }
00061 
00062 
00063 MeshRectangle* MeshRectangle::copy() const {
00064          MeshRectangle *returnvalue     = new MeshRectangle();
00065          
00066          returnvalue->start_        = this->start_; // apperently vector::operator=() takes a deep copy
00067          returnvalue->stop_         = this->stop_ ;
00068          returnvalue->multiplicity_ = this->multiplicity_;
00069          returnvalue->constDir_     = this->constDir_;
00070          return returnvalue;
00071 }
00072 
00073 int MeshRectangle::constDirection() const {
00074         return constDir_;
00075 }
00076 
00077 double MeshRectangle::constParameter() const {
00078         if(constDir_ > 2 || constDir_ < 0) return -1.0;
00079         return start_[constDir_];
00080 }
00081 
00082 int MeshRectangle::nKnotsIn(Basisfunction *basis) const {
00083         int hits = 0;
00084         int d    = constDir_;
00085         for(int i=0; i<=basis->getOrder(d); i++)
00086                 if( fabs((*basis)[d][i] - start_[d]) < DOUBLE_TOL )
00087                         hits++;
00088         return hits;
00089 }
00090 
00091 bool MeshRectangle::equals(const MeshRectangle *rect) const {
00092         for(int i=0; i<3; i++) {
00093                 if(start_[i] != rect->start_[i])     return false;
00094                 if(stop_[i]  != rect->stop_[i])      return false;
00095         }
00096         if(multiplicity_ != rect->multiplicity_) return false;
00097 
00098         return true;
00099 }
00100 
00102 bool MeshRectangle::overlaps(const MeshRectangle *rect) const {
00103         int c1 = constDir_; // constant index
00104         int v1 = (c1+1)%3;  // first variable index
00105         int v2 = (c1+2)%3;  // second variable index
00106 
00107         if(constDir_ == rect->constDir_) {
00108                 if(start_[c1] == rect->start_[c1] &&
00109                    start_[v1] <= rect->stop_[v1] && stop_[v1] >= rect->start_[v1] &&
00110                    start_[v2] <= rect->stop_[v2] && stop_[v2] >= rect->start_[v2])
00111                         return true;
00112         }
00113 
00114         return false;
00115 }
00116 
00118 bool MeshRectangle::contains(const MeshRectangle *rect) const {
00119         int c1 = constDir_; // constant index
00120         int v1 = (c1+1)%3;  // first variable index
00121         int v2 = (c1+2)%3;  // second variable index
00122 
00123         if(constDir_ == rect->constDir_) {
00124                 if(start_[c1] == rect->start_[c1] &&
00125                    start_[v1] <= rect->start_[v1] && stop_[v1] >= rect->stop_[v1] &&
00126                    start_[v2] <= rect->start_[v2] && stop_[v2] >= rect->stop_[v2] )
00127                         return true;
00128         }
00129 
00130         return false;
00131 }
00132 
00133 int MeshRectangle::makeOverlappingRects(std::vector<MeshRectangle*> &newGuys, int meshIndex, bool allowSplits) {
00134         int c1  = constDir_; // constant index
00135         int v1  = (c1+1)%3;  // first variable index
00136         int v2  = (c1+2)%3;  // second variable index
00137         int v[] = {v1, v2};  // index way of referencing these
00138         bool addThisToNewGuys = false;
00139         MeshRectangle *rect = newGuys.at(meshIndex);
00140         if( ! this->overlaps(rect) )
00141                 return 0;
00142         if( this->contains(rect) ) {
00143                 newGuys.erase(newGuys.begin() + meshIndex);
00144                 // std::cout << "Deleted: " << *rect << std::endl;
00145                 // std::cout << "  contained in : " << *this << std::endl;
00146                 delete rect;
00147                 return 1;
00148         }
00149         if( rect->contains(this) ) {
00150                 // std::cout << "Deleted: " << *this << std::endl;
00151                 // std::cout << "  contained in : " << *rect << std::endl;
00152                 return 2;
00153         }
00154 
00155 
00156         // for both variable indices... i=0 corresponds to checking everything left-right (direction v1),
00157         // i=1 is for checking everything up-down (direction v2). 
00158         for(int i=0; i<2; i++) {
00159                 int j = 1-i; // 0 or 1
00160                 /*
00161                  *    MERGE RECTANGLES TO ONE
00162                  *  +-------+------+
00163                  *  |       |      |
00164                  *  |       |      |
00165                  *  |       |      |
00166                  *  +-------+------+
00167                  */
00168                 // if the two mesh rectangles perfectly line up, keep only one of them
00169                 if((stop_[v[i]]  == rect->start_[v[i]] ||
00170                     start_[v[i]] == rect->stop_[v[i]]  )) {
00171 
00172                         if(start_[v[j]] == rect->start_[v[j]]  && 
00173                            stop_[v[j]]  == rect->stop_[v[j]]  ) {
00174                                 double min = std::min(start_[v[i]], rect->start_[v[i]]);
00175                                 double max = std::max(stop_[v[i]],  rect->stop_[v[i]] );
00176                                 // std::cout << "Deleted: " << *rect << std::endl;
00177                                 // std::cout << "  merged with  : " << *this << std::endl;
00178                                 newGuys.erase(newGuys.begin() + meshIndex);
00179                                 delete rect;
00180                                 start_[v[i]] = min;
00181                                 stop_[v[i]]  = max;
00182                                 if(addUniqueRect(newGuys, this))
00183                                         return 4;
00184                                 else
00185                                         return 5;
00186                         /*
00187                         *    ADD ELONGATED RECT
00188                         * y3 +-------+                
00189                         * y2 |       +------+         +-------+------+
00190                         *    |       |      |   =>    |              | 
00191                         * y1 +-------+      |         +-------+------+
00192                         *            |      |                new one
00193                         * y0         +------+               
00194                         *   x0      x1     x3
00195                         */
00196 
00197                         } else if((start_[v[j]] < rect->start_[v[j]]   && 
00198                                    stop_[v[j]]  < rect->stop_[v[j]] ) ||
00199                                   (start_[v[j]] > rect->start_[v[j]]   && 
00200                                   stop_[v[j]]  > rect->stop_[v[j]]  )) {
00201                                 double x0 = std::min(rect->start_[v[i]], start_[v[i]]);
00202                                 double x3 = std::max(rect->stop_[v[i]],  stop_[v[i]] );
00203                                 double y1 = std::max(rect->start_[v[j]], start_[v[j]]);
00204                                 double y2 = std::min(rect->stop_[v[j]],  stop_[v[j]] );
00205                                 double start[3];
00206                                 double stop[3];
00207                                 start[c1]   = start_[c1];
00208                                 stop[c1]    = stop_[c1];
00209                                 start[v[i]] = x0;
00210                                 stop[v[i]]  = x3;
00211                                 start[v[j]] = y1;
00212                                 stop[v[j]]  = y2;
00213                                 if(allowSplits) {
00214                                         MeshRectangle *m1 = new MeshRectangle(start, stop);
00215                                         if(!addUniqueRect(newGuys, m1))
00216                                                 delete m1;
00217                                 }
00218                         }
00219                 }
00220                 /*
00221                  *     EXPAND 'RECT' (RIGHT ONE)
00222                  *  +-------+
00223                  *  |    +--+------+
00224                  *  |    |  |      |
00225                  *  |    +--+------+
00226                  *  +-------+
00227                  */
00228                 if(start_[v[i]] <  rect->start_[v[i]] &&
00229                    start_[v[j]] <= rect->start_[v[j]] &&
00230                    stop_[v[j]]  >= rect->stop_[v[j]]) { // expand the support of rect DOWN in v[i]
00231                         rect->start_[v[i]] = start_[v[i]];
00232                 }
00233                 if(stop_[v[i]]  >  rect->stop_[v[i]]  &&
00234                    start_[v[j]] <= rect->start_[v[j]] &&
00235                    stop_[v[j]]  >= rect->stop_[v[j]]) { // expand the support of rect UP in v[i]
00236                         rect->stop_[v[i]] = stop_[v[i]];
00237                 }
00238                 /*
00239                  *     EXPAND 'THIS' (LEFT ONE)
00240                  *        +------+
00241                  *  +-----+--+   |
00242                  *  |     |  |   |
00243                  *  +-----+--+   |
00244                  *        +------+
00245                  */
00246                 if(rect->start_[v[i]] <  start_[v[i]] &&
00247                    rect->start_[v[j]] <= start_[v[j]] &&
00248                    rect->stop_[v[j]]  >= stop_[v[j]]) {
00249                         start_[v[i]] = rect->start_[v[i]];
00250                         addThisToNewGuys = true;
00251                 }
00252                 if(rect->stop_[v[i]]  >  stop_[v[i]]  &&
00253                    rect->start_[v[j]] <= start_[v[j]] &&
00254                    rect->stop_[v[j]]  >= stop_[v[j]]) {
00255                         stop_[v[i]] = rect->stop_[v[i]];
00256                         addThisToNewGuys = true;
00257                 }
00258         }
00259 
00260 
00261         /*
00262          *     DO NOTHING
00263          *          +--+               +--+          +----+      
00264          *          |  |               |  |          |    |
00265          *     +----+--+----+     +----+--+     +----+----+
00266          *     |    |  |    |     |    |  |     |    |    |
00267          *     +----+--+----+     +----+--+     |    |    |
00268          *          |  |               |  |     +----+----+ 
00269          *          +--+               +--+     
00270          *        note that this is after fixes above
00271          */
00272         if((rect->stop_[v1]  >=       stop_[v1]   &&
00273             rect->start_[v1] <=       start_[v1]  &&
00274                   stop_[v2]  >= rect->stop_[v2]   &&
00275                   start_[v2] <= rect->start_[v2]) ||
00276            (      stop_[v1]  >= rect->stop_[v1]   &&
00277                   start_[v1] <= rect->start_[v1]  &&
00278             rect->stop_[v2]  >=       stop_[v2]   &&
00279             rect->start_[v2] <=       start_[v2])) {
00280                 // ...
00281                 ;
00282         /*
00283          *         MAKE TWO NEW ONES
00284          *  y3  +-------+                        +--+
00285          *      |       |                        |  |
00286          *  y2  |    +--+----+           y2 +----+--+----+
00287          *      |    |  |    |      =>      |    |  |    | 
00288          *  y1  +----+--+    |           y1 +----+--+----+
00289          *           |       |                   |  |
00290          *  y0       +-------+                   +--+
00291          *     x0   x1 x2   x3                  x1  x2
00292          *                                    these are the two new ones
00293          */
00294         } else {
00295                 double x0 = std::min(rect->start_[v1], start_[v1]);
00296                 double x1 = std::max(rect->start_[v1], start_[v1]);
00297                 double x2 = std::min(rect->stop_[v1],  stop_[v1] );
00298                 double x3 = std::max(rect->stop_[v1],  stop_[v1] );
00299                 double y0 = std::min(rect->start_[v2], start_[v2]);
00300                 double y1 = std::max(rect->start_[v2], start_[v2]);
00301                 double y2 = std::min(rect->stop_[v2],  stop_[v2] );
00302                 double y3 = std::max(rect->stop_[v2],  stop_[v2] );
00303                 double start[3];
00304                 double stop[3];
00305                 start[c1] = start_[c1];
00306                 stop[c1]  = stop_[c1];
00307 
00308                 if(allowSplits) {
00309                         start[v1] = x0;
00310                         stop[v1]  = x3;
00311                         start[v2] = y1;
00312                         stop[v2]  = y2;
00313                         MeshRectangle *m1 = new MeshRectangle(start, stop);
00314 
00315                         start[v1] = x1;
00316                         stop[v1]  = x2;
00317                         start[v2] = y0;
00318                         stop[v2]  = y3;
00319                         MeshRectangle *m2 = new MeshRectangle(start, stop);
00320 
00321                         // std::cout << "Added: " << *m1 << std::endl;
00322                         // std::cout << "Added: " << *m2 << std::endl;
00323                         // std::cout << "  overlaps from  : " << *rect << std::endl;
00324                         // std::cout << "  overlaps from  : " << *this << std::endl;
00325 
00326                         if(!addUniqueRect(newGuys, m1))
00327                                 delete m1;
00328                         if(!addUniqueRect(newGuys, m2))
00329                                 delete m2;
00330                 }
00331 
00332         }
00333         if(addThisToNewGuys) {
00334                 // std::cout << "Moved: " << *this << std::endl;
00335                 if(addUniqueRect(newGuys, this)) 
00336                         return 3;
00337                 else
00338                         return 6;
00339         }
00340         return 0;
00341 }
00342 
00343 bool MeshRectangle::splits(Element *el) const {
00344         if(constDir_ == -1) {
00345                 std::cerr << "MeshRectangle::splits(Element*) - constDir_ not set on meshrectangle\n";
00346                 exit(3421);
00347         }
00348         
00349         int c1 = constDir_; // constant index
00350         int v1 = (c1+1)%3;  // first variable index
00351         int v2 = (c1+2)%3;  // second variable index
00352 
00353         if(el->getParmin(c1) < start_[c1]  && start_[c1] < el->getParmax(c1) && 
00354            start_[v1] <= el->getParmin(v1) && el->getParmax(v1) <= stop_[v1] && 
00355            start_[v2] <= el->getParmin(v2) && el->getParmax(v2) <= stop_[v2])
00356                 return true;
00357 
00358         return false;
00359 }
00360 
00361 bool MeshRectangle::splits(Basisfunction *basis) const {
00362         if(constDir_ == -1) {
00363                 std::cerr << "MeshRectangle::splits(Basisfunction*) - constDir_ not set on meshrectangle\n";
00364                 exit(3421);
00365         }
00366         
00367         int c1 = constDir_; // constant index
00368         int v1 = (c1+1)%3;  // first variable index
00369         int v2 = (c1+2)%3;  // second variable index
00370 
00371         if(basis->getParmin(c1) < start_[c1]  && start_[c1] < basis->getParmax(c1) && 
00372            start_[v1] <= basis->getParmin(v1) && basis->getParmax(v1) <= stop_[v1] && 
00373            start_[v2] <= basis->getParmin(v2) && basis->getParmax(v2) <= stop_[v2])
00374                 return true;
00375 
00376         return false;
00377 }
00378 
00379 bool MeshRectangle::operator==(const MeshRectangle &other) const {
00380         return this->equals(&other);
00381 }
00382 
00383 // convenience macro for reading formated input
00384 #define ASSERT_NEXT_CHAR(c) {ws(is); nextChar = is.get(); if(nextChar!=c) { std::cerr << "Error parsing meshrectangle\n"; std::cout << is; exit(325); } ws(is); }
00385 void MeshRectangle::read(std::istream &is) {
00386         char nextChar;
00387         ws(is);
00388 
00389         ASSERT_NEXT_CHAR('[');
00390         is >> start_[0];
00391         ASSERT_NEXT_CHAR(',');
00392         is >> stop_[0];
00393         ASSERT_NEXT_CHAR(']');
00394         ASSERT_NEXT_CHAR('x');
00395 
00396         ASSERT_NEXT_CHAR('[');
00397         is >> start_[1];
00398         ASSERT_NEXT_CHAR(',');
00399         is >> stop_[1];
00400         ASSERT_NEXT_CHAR(']');
00401         ASSERT_NEXT_CHAR('x');
00402 
00403         ASSERT_NEXT_CHAR('[');
00404         is >> start_[2];
00405         ASSERT_NEXT_CHAR(',');
00406         is >> stop_[2];
00407         ASSERT_NEXT_CHAR(']');
00408 
00409         ASSERT_NEXT_CHAR('(');
00410         is >> multiplicity_;
00411         ASSERT_NEXT_CHAR(')');
00412 
00413 
00414         for(int i=0; i<3; i++)
00415                 if(start_[i] == stop_[i])
00416                         constDir_ = i;
00417         if(constDir_ == -1)
00418                 std::cerr << "Error creating MeshRectangle: Not parallel to the parametric axis\n";
00419 }
00420 #undef ASSERT_NEXT_CHAR
00421 
00422 void MeshRectangle::write(std::ostream &os) const {
00423         os << "[" << start_[0] << ", " << stop_[0] << "] x ";
00424         os << "[" << start_[1] << ", " << stop_[1] << "] x ";
00425         os << "[" << start_[2] << ", " << stop_[2] << "] ";
00426         os << "(" << multiplicity_ << ")";
00427 }
00428 
00429 #undef DOUBLE_TOL
00430 
00431 } // end namespace LR
00432 
00433