Mesh.cc 28.5 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209
#include "AdaptStationary.h"
#include "AdaptInstationary.h"
#include "FiniteElemSpace.h"
#include "ElementData.h"
#include "MacroElement.h"
#include "MacroReader.h"
#include "Mesh.h"
#include "Traverse.h"
#include "Parameters.h"
#include "FixVec.h"
#include "DOFVector.h"
#include "CoarseningManager.h"
#include "DOFIterator.h"
#include "VertexVector.h"
#include "MacroWriter.h"
#include "PeriodicMap.h"
#include "Projection.h"
#include <algorithm>
#include <set>
#include <map>

#include "time.h"

namespace AMDiS {

#define TIME_USED(f,s) ((double)((s)-(f))/(double)CLOCKS_PER_SEC)

  //**************************************************************************
  //  flags, which information should be present in the elInfo structure     
  //**************************************************************************

  const Flag Mesh::FILL_NOTHING    = 0X00L;
  const Flag Mesh::FILL_COORDS     = 0X01L;
  const Flag Mesh::FILL_BOUND      = 0X02L;
  const Flag Mesh::FILL_NEIGH      = 0X04L;
  const Flag Mesh::FILL_OPP_COORDS = 0X08L;
  const Flag Mesh::FILL_ORIENTATION= 0X10L;
  const Flag Mesh::FILL_DET        = 0X20L;
  const Flag Mesh::FILL_GRD_LAMBDA = 0X40L;
  const Flag Mesh::FILL_ADD_ALL    = 0X80L;


  const Flag Mesh::FILL_ANY_1D= (0X01L|0X02L|0X04L|0X08L|0x20L|0X40L|0X80L);
  const Flag Mesh::FILL_ANY_2D= (0X01L|0X02L|0X04L|0X08L|0x20L|0X40L|0X80L);
  const Flag Mesh::FILL_ANY_3D= (0X01L|0X02L|0X04L|0X08L|0X10L|0x20L|0X40L|0X80L);

  //**************************************************************************
  //  flags for Mesh traversal                                                
  //**************************************************************************

  const Flag Mesh::CALL_EVERY_EL_PREORDER  = 0X0100L;
  const Flag Mesh::CALL_EVERY_EL_INORDER   = 0X0200L;
  const Flag Mesh::CALL_EVERY_EL_POSTORDER = 0X0400L;
  const Flag Mesh::CALL_LEAF_EL            = 0X0800L;
  const Flag Mesh::CALL_LEAF_EL_LEVEL      = 0X1000L;
  const Flag Mesh::CALL_EL_LEVEL           = 0X2000L;
  const Flag Mesh::CALL_MG_LEVEL           = 0X4000L ; // used in mg methods 


  // const Flag Mesh::USE_PARAMETRIC          = 0X8000L ; // used in mg methods 

  // ::std::list<Mesh*> Mesh::meshes;
  DOFAdmin* Mesh::compressAdmin = NULL;
  Mesh* Mesh::traversePtr = NULL;
  int Mesh::iadmin = 0;
  ::std::vector<DegreeOfFreedom> Mesh::dof_used;
  const int Mesh::MAX_DOF=100;
  ::std::map<DegreeOfFreedom, DegreeOfFreedom*> Mesh::serializedDOFs;

  struct delmem { 
    DegreeOfFreedom* ptr;
    int              len;
  };


  Mesh::Mesh(const ::std::string& aName, int dimension) 
    : name(aName), 
      dim(dimension), 
      nVertices(0),
      nEdges(0),
      nLeaves(0), 
      nElements(0),
      parametric(NULL), 
      preserveCoarseDOFs(false),
      nDOFEl(0),
      nDOF(dimension, DEFAULT_VALUE, 0),
      nNodeEl(0),
      node(dimension, DEFAULT_VALUE, 0),
      elementPrototype(NULL),
      elementDataPrototype(NULL),
      elementIndex(-1),
      initialized(false),
      final_lambda(dimension, DEFAULT_VALUE, 0.0)
  {

    FUNCNAME("Mesh::Mesh");

    // set default element prototype
    switch(dim) {
    case 1:
      elementPrototype = NEW Line(this);
      break;
    case 2:
      elementPrototype = NEW Triangle(this);
      break;
    case 3:
      elementPrototype = NEW Tetrahedron(this);
      break;
    default:
      ERROR_EXIT("invalid dimension\n");
    }

    elementPrototype->setIndex(-1);

    elementIndex=0;
  };

  Mesh::~Mesh()
  {
  };

  void Mesh::addMacroElement(MacroElement* m) {
    macroElements.push_back(m); 
    m->setIndex(macroElements.size());
  };




  int Mesh::traverse(int level, Flag flag, 
		     int (*el_fct)(ElInfo*))
  {
    FUNCNAME("Mesh::traverse()");
    ::std::deque<MacroElement*>::iterator mel;
    ElInfo* elinfo = createNewElInfo();
    Traverse tinfo(this, flag, level, el_fct);
    int sum = 0;
  
    elinfo->setFillFlag(flag);
  
    if (flag.isSet(Mesh::CALL_LEAF_EL_LEVEL) || 
	flag.isSet(Mesh::CALL_EL_LEVEL)      || 
	flag.isSet(Mesh::CALL_MG_LEVEL)) {
      TEST(level >= 0)("invalid level: %d\n", level);
    }
  
    for (mel = macroElements.begin(); mel != macroElements.end(); mel++) {
      elinfo->fillMacroInfo(*mel);
      sum += tinfo.recursive(elinfo);
    }

    DELETE elinfo;
    
    return (flag.isSet(Mesh::FILL_ADD_ALL)) ? sum : 0;
  }



  void Mesh::addDOFAdmin(DOFAdmin *localAdmin)
  {    
    FUNCNAME("Mesh::addDOFAdmin()");

    int i, j, d, n;
    ::std::vector<DOFAdmin*>::iterator dai;

    localAdmin->setMesh(this);
    n = admin.size();

    dai=::std::find(admin.begin(),admin.end(),localAdmin);
    if (dai!= admin.end()) {
      ERROR("admin %s is already associated to mesh %s\n",
	    localAdmin->getName().c_str(), this->getName().c_str());
    }

    // ===== adding dofs to already existing elements ============================ 

    if (initialized) {
      static bool pnd_1d_0[2] = {true, true};
      static bool pnd_1d_1[1] = {false};
      static bool pnd_2d_0[3] = {true, true, true};
      static bool pnd_2d_1[3] = {true, true, false};
      static bool pnd_2d_2[1] = {false};
      static bool pnd_3d_0[4] = {true, true, true, true};
      static bool pnd_3d_1[6] = {false, true, true, true, true, true};
      static bool pnd_3d_2[4] = {true, true, false, false};
      static bool pnd_3d_3[1] = {false};
      static bool *pnd_1d[2] = {pnd_1d_0, pnd_1d_1};
      static bool *pnd_2d[3] = {pnd_2d_0, pnd_2d_1, pnd_2d_2};
      static bool *pnd_3d[4] = {pnd_3d_0, pnd_3d_1, pnd_3d_2, pnd_3d_3};
      static bool **parentNeedsDOF[4] = {NULL, pnd_1d, pnd_2d, pnd_3d};

     
      ::std::list<struct delmem> delList;
      ::std::map< ::std::set<DegreeOfFreedom>, DegreeOfFreedom*> dofPtrMap;
      const DOFAdmin *vertexAdmin = getVertexAdmin();
      int vertexAdminPreDOFs = vertexAdmin->getNumberOfPreDOFs(VERTEX);

      // finding necessary node number for new admin

      int newNNode=0;
      GeoIndex geoIndex;

      for(d = 0; d < dim+1; d++) {
	geoIndex = INDEX_OF_DIM(d, dim);
      
	if (localAdmin->getNumberOfDOFs(geoIndex)>0||nDOF[geoIndex]>0)
	  newNNode+=getGeo(geoIndex);
      };

210 211
      bool extendNodes = (newNNode>nNodeEl);
      int oldNNodes = nNodeEl;
212

213
      nNodeEl = newNNode;
214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246

      TraverseStack stack;
      ElInfo *elInfo = NULL;
    
      WARNING("You are using untested code (adding dofs to existing mesh). Please contact\nsoftware administrator if any errors occur in this context.\n");

      elInfo = stack.traverseFirst(this, -1, CALL_EVERY_EL_PREORDER);
      while(elInfo) {
	Element *element = elInfo->getElement();
	DegreeOfFreedom *newDOF, **oldDOF, **dof = 
	  const_cast<DegreeOfFreedom**>(element->getDOF());

	int index = 0;

	if (extendNodes) {
	  oldDOF=dof;
	  element->setDOFPtrs();
	  dof=const_cast<DegreeOfFreedom**>(element->getDOF());
	  int index=0,oldIndex=0;
	  for(d = 0; d < dim+1; d++) {
	    geoIndex = INDEX_OF_DIM(d, dim);
	    if (nDOF[geoIndex]>0) {
	      for(i=0;i<getGeo(geoIndex);++i) 
		dof[index++]=oldDOF[oldIndex++];
	    }
	    else {
	      if (localAdmin->getNumberOfDOFs(geoIndex)>0) 
		index+=getGeo(geoIndex);
	    }
	  }
	
	  FREE_MEMORY(oldDOF, DegreeOfFreedom*, oldNNodes);

247
	  TEST_EXIT_DBG(index == nNodeEl)("ERROR: Number of entered nodes %f != number of nodes %f\n",index,nNodeEl);
248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334

	}


	index=0;

	// allocate new memory at elements
	for(d = 0; d < dim+1; d++) {
	  geoIndex = INDEX_OF_DIM(d, dim);
      
	  int numberOfDOFs = localAdmin->getNumberOfDOFs(geoIndex);
	  int numberOfPreDOFs = nDOF[geoIndex];

	  if (numberOfDOFs>0||numberOfPreDOFs>0) {

	    // for all vertices/edges/...
	    for(i = 0; i < getGeo(geoIndex); i++, index++) {
	      ::std::set<DegreeOfFreedom> dofSet;
	      for(j = 0; j < d+1; j++) {
		dofSet.insert(dof[element->getVertexOfPosition(geoIndex, i, j)][vertexAdminPreDOFs]);
	      }
	    
	      if(element->isLeaf() || parentNeedsDOF[dim][d][i]) {
		if(dofPtrMap[dofSet] == NULL) {
		  if(localAdmin->getNumberOfDOFs(geoIndex)) {
		    newDOF = GET_MEMORY(DegreeOfFreedom, numberOfPreDOFs + numberOfDOFs);
		    // copy old dofs to new memory and free old memory
		    if(dof[index]) {
		      for(j = 0; j < numberOfPreDOFs; j++) {
			newDOF[j] = dof[index][j];
		      }
		      //	  FREE_MEMORY(dof[index], DegreeOfFreedom, numberOfPreDOFs);
		      // Do not free memory. The information has to be used to identify the part in other elements.
		      // The memory is only marked for freeing.
		      struct delmem fm;
		      fm.ptr=dof[index];
		      fm.len=numberOfPreDOFs;
		      delList.push_back(fm);
		    }
		    for(j = 0; j < numberOfDOFs; j++) {
		      newDOF[numberOfPreDOFs + j] = localAdmin->getDOFIndex();
		    }
		    dof[index] = newDOF;
		  }
		  dofPtrMap[dofSet] = dof[index];
		} else {
		  dof[index] = dofPtrMap[dofSet];
		}
	      }
	    }
	  }
	}
	elInfo = stack.traverseNext(elInfo);
      }
  
      // now free the old dof memory:

      ::std::list<struct delmem>::iterator it=delList.begin();
    
      while(it!=delList.end()) {
	FREE_MEMORY((*it).ptr, DegreeOfFreedom, (*it).len);
	it++;
      }

      delList.clear();

    }
    // ============================================================================

    admin.push_back(localAdmin);

    nDOFEl = 0;

    localAdmin->setNumberOfPreDOFs(VERTEX,nDOF[VERTEX]);
    nDOF[VERTEX]  += localAdmin->getNumberOfDOFs(VERTEX);
    nDOFEl += getGeo(VERTEX) * nDOF[VERTEX];

    if(dim > 1) {
      localAdmin->setNumberOfPreDOFs(EDGE,nDOF[EDGE]);
      nDOF[EDGE]    += localAdmin->getNumberOfDOFs(EDGE);
      nDOFEl += getGeo(EDGE) * nDOF[EDGE];
    }

    localAdmin->setNumberOfPreDOFs(CENTER,nDOF[CENTER]);
    nDOF[CENTER]  += localAdmin->getNumberOfDOFs(CENTER);
    nDOFEl += nDOF[CENTER];

335
    TEST_EXIT_DBG(nDOF[VERTEX] > 0)("no vertex dofs\n");
336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365

    node[VERTEX]  = 0;
    nNodeEl     = getGeo(VERTEX);

    if(dim > 1) {
      node[EDGE]    = nNodeEl;
      if (nDOF[EDGE] > 0) nNodeEl += getGeo(EDGE);
    }

    if (3==dim){
      localAdmin->setNumberOfPreDOFs(FACE,nDOF[FACE]);
      nDOF[FACE]  += localAdmin->getNumberOfDOFs(FACE);
      nDOFEl     += getGeo(FACE) * nDOF[FACE];
      node[FACE]    = nNodeEl;
      if (nDOF[FACE] > 0) nNodeEl +=  getGeo(FACE);
    }

    node[CENTER]    = nNodeEl;
    if (nDOF[CENTER] > 0) nNodeEl += 1;

    return;
  }


  /****************************************************************************/
  /*  dofCompress: remove holes in dof vectors                                */
  /****************************************************************************/

  void Mesh::dofCompress()
  {
366 367 368
    FUNCNAME("Mesh::dofCompress()");
    int size;
    Flag fill_flag;
369

370
    for (iadmin = 0; iadmin < static_cast<int>(admin.size()); iadmin++) {
371 372 373
      compressAdmin = admin[iadmin];

      TEST_EXIT_DBG(compressAdmin)("no admin[%d] in mesh\n", iadmin);
374 375 376
      
      if ((size = compressAdmin->getSize()) < 1) 
	continue;
Thomas Witkowski's avatar
Thomas Witkowski committed
377

378 379
      if (compressAdmin->getUsedDOFs() < 1)    
	continue;
Thomas Witkowski's avatar
Thomas Witkowski committed
380

381 382
      if (compressAdmin->getHoleCount() < 1)    
	continue;
Thomas Witkowski's avatar
Thomas Witkowski committed
383
  
384 385 386 387 388 389 390 391
      newDOF.resize(size);
      
      compressAdmin->compress(newDOF);
      
      if (preserveCoarseDOFs) {
	fill_flag = Mesh::CALL_EVERY_EL_PREORDER | Mesh::FILL_NOTHING;
      } else {
	fill_flag = Mesh::CALL_LEAF_EL | Mesh::FILL_NOTHING;
392
      }
393 394 395 396 397 398
      
      traverse( -1, fill_flag, newDOFFct1);
      traverse( -1, fill_flag, newDOFFct2);
      
      newDOF.resize(0);
    }   
399 400 401 402 403
  }


  DegreeOfFreedom *Mesh::getDOF(GeoIndex position)
  {
404
    FUNCNAME("Mesh::getDOF()");
405

406
    TEST_EXIT_DBG(position >= CENTER && position <= FACE)
407
      ("unknown position %d\n", position);
408

409 410 411
    int ndof = getNumberOfDOFs(position);
    if (ndof <= 0) 
      return(NULL);
412

413
    DegreeOfFreedom *dof = GET_MEMORY(DegreeOfFreedom, ndof);
414

415 416
    for (int i = 0; i < getNumberOfDOFAdmin(); i++) {
      const DOFAdmin *localAdmin = &getDOFAdmin(i);
417
      TEST_EXIT_DBG(localAdmin)("no admin[%d]\n", i);
418 419 420 421
      
      int n  = localAdmin->getNumberOfDOFs(position);
      int n0 = localAdmin->getNumberOfPreDOFs(position);
      
422
      TEST_EXIT_DBG(n + n0 <= ndof)("n=%d, n0=%d too large: ndof=%d\n", n, n0, ndof);
423 424 425
      
      for (int j = 0; j < n; j++) {
	dof[n0 + j] = const_cast<DOFAdmin*>(localAdmin)->getDOFIndex();
426
      }
427
    }
428 429 430 431 432 433 434
  
    return(dof);
  }


  DegreeOfFreedom **Mesh::createDOFPtrs()
  {
435
    FUNCNAME("Mesh::createDOFPtrs()");
436 437 438 439

    if (nNodeEl <= 0)
      return(NULL);

440 441
    DegreeOfFreedom **ptrs = GET_MEMORY(DegreeOfFreedom*, nNodeEl);
    for (int i = 0; i < nNodeEl; i++)
442 443 444 445 446 447 448
      ptrs[i] = NULL;

    return(ptrs);
  }

  void Mesh::freeDOFPtrs(DegreeOfFreedom **ptrs)
  {
449
    FUNCNAME("Mesh::freeDOFPtrs()");
450

451
    TEST_EXIT_DBG(ptrs)("ptrs=NULL\n");
452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523

    if (nNodeEl <= 0)
      return;
  
    FREE_MEMORY(ptrs, DegreeOfFreedom*, nNodeEl);
  }


  const DOFAdmin *Mesh::createDOFAdmin(const ::std::string& lname,DimVec<int> lnDOF)
  {
    FUNCNAME("Mesh::createDOFAdmin");

    DOFAdmin         *localAdmin;
    int              i;

    localAdmin=NEW DOFAdmin(this,lname);

    for (i = 0; i < dim+1; i++)
      localAdmin->setNumberOfDOFs(i,lnDOF[i]);

    addDOFAdmin(localAdmin);

    return(localAdmin);
  }





  // int Mesh::macroType(const ::std::string& filename, const ::std::string& type)
  // {
  //   const char *fn, *t;

  //   if (3==dim) return 0;
  
  //   if (filename.size() <= type.size())
  //     return(false);

  //   fn = filename.data();
  //   while (*fn) fn++;
  //   t = type.data();
  //   while (*t) t++;

  //   while (t != type  &&  *t == *fn) t--;
  
  //   return(t == type);
  // }

  const DOFAdmin* Mesh::getVertexAdmin() const
  {
    int       i;
    const DOFAdmin *localAdmin = NULL;

    for (i = 0; i < static_cast<int>(admin.size()); i++)
      {
	if (admin[i]->getNumberOfDOFs(VERTEX))
	  {
	    if (!localAdmin)  
	      localAdmin = admin[i];
	    else if (admin[i]->getSize() < localAdmin->getSize())
	      localAdmin = admin[i];
	  }
      }
    return(localAdmin);
  }

  void Mesh::freeDOF(DegreeOfFreedom* dof, GeoIndex position)
  {
    FUNCNAME("Mesh::freeDOF");
    DOFAdmin *localAdmin;
    int     i, j, n, n0, ndof;

524
    TEST_EXIT_DBG(position >= CENTER && position <= FACE)
525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544
      ("unknown position %d\n",position);

    ndof = nDOF[position];
    if (ndof) 
      {
	if (!dof)
	  {
	    MSG("dof = NULL, but ndof=%d\n", ndof);
	    return;
	  }
      }
    else
      {
	if (dof)
	  {
	    MSG("dof != NULL, but ndof=0\n");
	  }
	return;
      }

545
    TEST_EXIT_DBG(ndof <= MAX_DOF)
546 547 548 549 550 551 552 553 554
      ("ndof too big: ndof=%d, MAX_DOF=%d\n",ndof,MAX_DOF);

    for (i = 0; i < static_cast<int>(admin.size()); i++)
      {
	localAdmin = admin[i];

	n  = localAdmin->getNumberOfDOFs(position);
	n0 = localAdmin->getNumberOfPreDOFs(position);

555
	TEST_EXIT_DBG(n+n0 <= ndof)("n=%d, n0=%d too large: ndof=%d\n", n, n0, ndof);
556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576

	for (j = 0; j < n; j++)
	  {
	    localAdmin->freeDOFIndex(dof[n0+j]);
	  }
      }

    FREE_MEMORY(dof, DegreeOfFreedom, ndof);
    return;  
  }

  void Mesh::freeElement(Element* el)
  {
    freeDOFPtrs(const_cast<DegreeOfFreedom**>(el->getDOF()));
    DELETE el;
  }


  Element* Mesh::createNewElement(Element *parent)
  {
    FUNCNAME("Mesh::createNewElement()");
577 578

    TEST_EXIT_DBG(elementPrototype)("no element prototype\n");
579 580 581

    Element *el = parent ? parent->clone() : elementPrototype->clone();
  
582
    if (!parent && elementDataPrototype) {
583 584 585 586 587 588 589 590
      el->setElementData(elementDataPrototype->clone()); 
    } else {
      el->setElementData(NULL); // must be done in ElementData::refineElementData()
    }

    return el;
  }

591

592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627
  ElInfo* Mesh::createNewElInfo()
  {
    switch(dim) {
    case 1:
      return NEW ElInfo1d(this);
      break;
    case 2:
      return NEW ElInfo2d(this);
      break;
    case 3:
      return NEW ElInfo3d(this);
      break;
    default:
      ERROR_EXIT("invalid dim\n");
      return NULL;
    };
  }



  bool Mesh::findElInfoAtPoint(const WorldVector<double>& xy,
			       ElInfo *el_info,
			       DimVec<double>&    bary,
			       const MacroElement      *start_mel,
			       const WorldVector<double> *xy0,
			       double            *sp)
  {
    static const MacroElement *mel = NULL;
    DimVec<double> lambda(dim, NO_INIT);
    ElInfo *mel_info = NULL;

    mel_info = createNewElInfo();

    if (start_mel != NULL)
      mel = start_mel;
    else
628
      if ((mel == NULL) || (mel->getElement()->getMesh() != this))
629 630 631
	mel = *(macroElements.begin());

    mel_info->setFillFlag(Mesh::FILL_COORDS);
632
    g_xy = &xy;
633
    g_xy0 = xy0;
634
    g_sp = sp;
635 636 637

    mel_info->fillMacroInfo(mel);

638
    int k;
639 640 641 642 643 644 645 646 647 648
    while ((k = mel_info->worldToCoord(xy, &lambda)) >= 0) {
      if (mel->getNeighbour(k)) {
	mel = mel->getNeighbour(k);
	mel_info->fillMacroInfo(mel);
	continue;
      }
      break;
    }

    /* now, descend in tree to find leaf element at point */
649 650 651 652
    bool inside = findElementAtPointRecursive(mel_info, lambda, k, el_info);
    for (int i = 0; i <= dim; i++) {
      bary[i] = final_lambda[i];
    }
653 654 655 656 657 658 659
  
    DELETE mel_info;

    return(inside);
  }

  bool Mesh::findElementAtPoint(const WorldVector<double>&  xy,
660 661
				Element **elp, 
				DimVec<double>& bary,
662
				const MacroElement *start_mel,
663 664
				const WorldVector<double> *xy0,
				double *sp)
665
  {
666 667
    ElInfo *el_info = createNewElInfo();
    int val = findElInfoAtPoint(xy, el_info, bary, start_mel, xy0, sp);
668 669 670 671 672 673 674 675 676 677

    *elp = el_info->getElement();

    DELETE el_info;

    return(val);
  }



678
  bool Mesh::findElementAtPointRecursive(ElInfo *el_info,
679
					 const DimVec<double>& lambda,
680
					 int outside,
681 682
					 ElInfo* final_el_info)
  {
683
    FUNCNAME("Mesh::findElementAtPointRecursive()");
684 685
    Element *el = el_info->getElement();
    DimVec<double> c_lambda(dim, NO_INIT);
686 687
    int inside;
    int ichild, c_outside;
688 689 690 691

    if (el->isLeaf()) {
      *final_el_info = *el_info;
      if (outside < 0) {
692 693 694 695
	for (int i = 0; i <= dim; i++) {
	  final_lambda[i] = lambda[i];
	}

696
	return(true);
697 698 699 700 701 702 703 704 705 706 707 708
      }  else {  /* outside */
	if (g_xy0) { /* find boundary point of [xy0, xy] */
	  el_info->worldToCoord(*(g_xy0), &c_lambda);
	  double s = lambda[outside] / (lambda[outside] - c_lambda[outside]);
	  for (int i = 0; i <= dim; i++) {
	    final_lambda[i] = s * c_lambda[i] + (1.0-s) * lambda[i];
	  }
	  if (g_sp) {
	    *(g_sp) = s;
	  }
	  if (dim == 3) 
	    MSG("outside finest level on el %d: s=%.3e\n", el->getIndex(), s);
709

710
	  return(false);  /* ??? */
711
	}
712 713
	else return(false);
      }
714 715
    }

716
    ElInfo *c_el_info = createNewElInfo();
717

718
    if (dim == 1) {
719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739
      if (lambda[0] >= lambda[1]) {
	c_el_info->fillElInfo(0, el_info);
	if (outside >= 0) {
	  outside = el_info->worldToCoord(*(g_xy), &c_lambda);
	  if (outside >= 0) ERROR("point outside domain\n");
	} else {
	  c_lambda[0] = lambda[0] - lambda[1];
	  c_lambda[1] = 2.0 * lambda[1];
	}
      } else {
	c_el_info->fillElInfo(1, el_info);
	if (outside >= 0)  {
	  outside = el_info->worldToCoord(*(g_xy), &c_lambda);
	  if (outside >= 0) ERROR("point outside domain\n");
	} else {
	  c_lambda[1] = lambda[1] - lambda[0];
	  c_lambda[0] = 2.0 * lambda[0];
	}
      }
    } /* DIM == 1 */

740
    if (dim == 2) {
741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767
      if (lambda[0] >= lambda[1]) {
	c_el_info->fillElInfo(0, el_info);
	if (el->isNewCoordSet()) {
	  outside = c_el_info->worldToCoord(*(g_xy), &c_lambda);
	  if (outside >= 0) {
	    ERROR("outside curved boundary child 0\n");
	  }
	} else {
	  c_lambda[0] = lambda[2];
	  c_lambda[1] = lambda[0] - lambda[1];
	  c_lambda[2] = 2.0 * lambda[1];
	}
      } else {
	c_el_info->fillElInfo(1, el_info);
	if (el->isNewCoordSet()) {
	  outside = c_el_info->worldToCoord(*(g_xy), &c_lambda);
	  if (outside >= 0) {
	    ERROR("outside curved boundary child 1\n");
	  }
	} else {
	  c_lambda[0] = lambda[1] - lambda[0];
	  c_lambda[1] = lambda[2];
	  c_lambda[2] = 2.0 * lambda[0];
	}
      }
    } /* DIM == 2 */

768
    if (dim == 3) {
769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791
      if (el->isNewCoordSet()) {
	if (lambda[0] >= lambda[1])
	  ichild = 0;
	else
	  ichild = 1;
	c_el_info->fillElInfo(ichild, el_info);
	c_outside = c_el_info->worldToCoord(*(g_xy), &c_lambda);

	if (c_outside>=0) {  /* test is other child is better... */
	  DimVec<double> c_lambda2(dim, NO_INIT);
	  int c_outside2;
	  ElInfo *c_el_info2 = createNewElInfo();

	  c_el_info2->fillElInfo(1-ichild, el_info);
	  c_outside2 = c_el_info2->worldToCoord(*(g_xy), &c_lambda2);

	  MSG("new_coord CHILD %d: outside=%d, lambda=(%.2f %.2f %.2f %.2f)\n",
	      ichild, c_outside, c_lambda[0],c_lambda[1],c_lambda[2],c_lambda[3]);
	  MSG("new_coord CHILD %d: outside=%d, lambda=(%.2f %.2f %.2f %.2f)\n",
	      1-ichild, c_outside2, c_lambda2[0],c_lambda2[1],c_lambda2[2],
	      c_lambda2[3]);

	  if ((c_outside2 < 0) || (c_lambda2[c_outside2] > c_lambda[c_outside])) {
792 793 794
	    for (int i = 0; i <= dim; i++) {
	      c_lambda[i] = c_lambda2[i];
	    }
795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925
	    c_outside = c_outside2;
	    *c_el_info = *c_el_info2;
	    ichild = 1 - ichild;
	  }
	  DELETE c_el_info2;
	}
	outside = c_outside;
      } else {  /* no new_coord */
	if (lambda[0] >= lambda[1]) {
	  c_el_info->fillElInfo(0, el_info);
	  c_lambda[0] = lambda[0] - lambda[1];
	  c_lambda[1] = lambda[Tetrahedron::childVertex[(dynamic_cast<ElInfo3d*>(el_info))->
							getType()][0][1]];
	  c_lambda[2] = lambda[Tetrahedron::childVertex[(dynamic_cast<ElInfo3d*>(el_info))->
							getType()][0][2]];
	  c_lambda[3] = 2.0 * lambda[1];
	} else {
	  c_el_info->fillElInfo(1, el_info);
	  c_lambda[0] = lambda[1] - lambda[0];
	  c_lambda[1] = lambda[Tetrahedron::childVertex[(dynamic_cast<ElInfo3d*>(el_info))->
							getType()][1][1]];
	  c_lambda[2] = lambda[Tetrahedron::childVertex[(dynamic_cast<ElInfo3d*>(el_info))->
							getType()][1][2]];
	  c_lambda[3] = 2.0 * lambda[0];
	}
      }
    }  /* DIM == 3 */

    inside = findElementAtPointRecursive(c_el_info, c_lambda, outside, 
					 final_el_info);
    DELETE c_el_info;

    return(inside); 
  }


  void Mesh::setDiameter(const WorldVector<double>& w) { diam = w; }

  void Mesh::setDiameter(int i, double w) { diam[i] = w; }


  int Mesh::newDOFFct1(ElInfo* ei) {
    ei->getElement()->newDOFFct1(compressAdmin);
    return 0;
  }

  int Mesh::newDOFFct2(ElInfo* ei) {
    ei->getElement()->newDOFFct2(compressAdmin);
    return 0;
  }

  void Mesh::serialize(::std::ostream &out)
  {
    serializedDOFs.clear();

    // write name
    out << name << ::std::endl;

    // write dim
    out.write(reinterpret_cast<const char*>(&dim), sizeof(int));

    // write nVertices
    out.write(reinterpret_cast<const char*>(&nVertices), sizeof(int));

    // write nEdges
    out.write(reinterpret_cast<const char*>(&nEdges), sizeof(int));

    // write nLeaves
    out.write(reinterpret_cast<const char*>(&nLeaves), sizeof(int));

    // write nElements
    out.write(reinterpret_cast<const char*>(&nElements), sizeof(int));

    // write nFaces
    out.write(reinterpret_cast<const char*>(&nFaces), sizeof(int));

    // write maxEdgeNeigh
    out.write(reinterpret_cast<const char*>(&maxEdgeNeigh), sizeof(int));

    // write diam
    diam.serialize(out);

    // write preserveCoarseDOFs
    out.write(reinterpret_cast<const char*>(&preserveCoarseDOFs), sizeof(bool));

    // write nDOFEl
    out.write(reinterpret_cast<const char*>(&nDOFEl), sizeof(int));

    // write nDOF
    nDOF.serialize(out);

    // write nNodeEl
    out.write(reinterpret_cast<const char*>(&nNodeEl), sizeof(int));

    // write node
    node.serialize(out);

    // write admins
    int i, size = static_cast<int>(admin.size());
    out.write(reinterpret_cast<const char*>(&size), sizeof(int));
    for (i = 0; i < size; i++) {
      admin[i]->serialize(out);
    }

    // write macroElements
    size = static_cast<int>(macroElements.size());
    out.write(reinterpret_cast<const char*>(&size), sizeof(int));
    for (i = 0; i < size; i++) {
      macroElements[i]->serialize(out);
    }

    // write elementIndex
    out.write(reinterpret_cast<const char*>(&elementIndex), sizeof(int));

    // write initialized
    out.write(reinterpret_cast<const char*>(&initialized), sizeof(bool));

    serializedDOFs.clear();
  }

  void Mesh::deserialize(::std::istream &in)
  {
    serializedDOFs.clear();

    // read name
    in >> name;
    in.get();

    // read dim
    int oldVal = dim;
    in.read(reinterpret_cast<char*>(&dim), sizeof(int));
926
    TEST_EXIT_DBG((oldVal == 0) || (dim == oldVal))("invalid dimension\n");
927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954

    // read nVertices
    in.read(reinterpret_cast<char*>(&nVertices), sizeof(int));

    // read nEdges
    in.read(reinterpret_cast<char*>(&nEdges), sizeof(int));

    // read nLeaves
    in.read(reinterpret_cast<char*>(&nLeaves), sizeof(int));

    // read nElements
    in.read(reinterpret_cast<char*>(&nElements), sizeof(int));

    // read nFaces
    in.read(reinterpret_cast<char*>(&nFaces), sizeof(int));

    // read maxEdgeNeigh
    in.read(reinterpret_cast<char*>(&maxEdgeNeigh), sizeof(int));

    // diam
    diam.deserialize(in);

    // read preserveCoarseDOFs
    in.read(reinterpret_cast<char*>(&preserveCoarseDOFs), sizeof(bool));

    // read nDOFEl
    oldVal = nDOFEl;
    in.read(reinterpret_cast<char*>(&nDOFEl), sizeof(int));
955
    TEST_EXIT_DBG((oldVal == 0) || (nDOFEl == oldVal))("invalid nDOFEl\n");
956 957 958 959 960 961 962

    // read nDOF
    nDOF.deserialize(in);

    // read nNodeEl
    oldVal = nNodeEl;
    in.read(reinterpret_cast<char*>(&nNodeEl), sizeof(int));
963
    TEST_EXIT_DBG((oldVal == 0) || (nNodeEl == oldVal))("invalid nNodeEl\n");
964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094

    // read node
    node.deserialize(in);

    // read admins
    int i, size;
    in.read(reinterpret_cast<char*>(&size), sizeof(int));
    admin.resize(size, NULL);
    for (i = 0; i < size; i++) {
      if (!admin[i]) {
	admin[i] = NEW DOFAdmin(this);
      }
      admin[i]->deserialize(in);
    }

    // read macroElements
    in.read(reinterpret_cast<char*>(&size), sizeof(int));

    ::std::vector< ::std::vector<int> > neighbourIndices(size);

    for (i = 0; i < static_cast<int>(macroElements.size()); i++) {
      if (macroElements[i]) {
	DELETE macroElements[i];
      }
    }
    macroElements.resize(size);
    for(i = 0; i < size; i++) {
      macroElements[i] = NEW MacroElement(dim);
      macroElements[i]->writeNeighboursTo(&(neighbourIndices[i]));
      macroElements[i]->deserialize(in);
    }

    // read elementIndex
    in.read(reinterpret_cast<char*>(&elementIndex), sizeof(int));

    // read initialized
    in.read(reinterpret_cast<char*>(&initialized), sizeof(bool));

    // set neighbour pointer in macro elements
    int j, neighs = getGeo(NEIGH);
    for(i = 0; i < static_cast<int>(macroElements.size()); i++) {
      for(j = 0; j < neighs; j++) {
	int index = neighbourIndices[i][j];
	if(index != -1) {
	  macroElements[i]->setNeighbour(j, macroElements[index]);
	} else {
	  macroElements[i]->setNeighbour(j, NULL);
	}
      }
    }

    // set mesh pointer in elements
    TraverseStack stack;
    ElInfo *elInfo = stack.traverseFirst(this, -1, CALL_EVERY_EL_PREORDER);
    while(elInfo) {
      elInfo->getElement()->setMesh(this);
      elInfo = stack.traverseNext(elInfo);
    }

    serializedDOFs.clear();
  }

  void Mesh::initialize() 
  {
    ::std::string macroFilename("");
    ::std::string valueFilename("");
    ::std::string periodicFile("");
    int check = 1;

    GET_PARAMETER(0, name + "->macro file name",  &macroFilename);
    GET_PARAMETER(0, name + "->value file name",  &valueFilename);
    GET_PARAMETER(0, name + "->periodic file", &periodicFile);
    GET_PARAMETER(0, name + "->check", "%d", &check);
    GET_PARAMETER(0, name + "->preserve coarse dofs", "%d", &preserveCoarseDOFs);

    if (macroFilename.length()) {
      macroFileInfo_ = MacroReader::readMacro(macroFilename.c_str(), 
					      this,
					      periodicFile == "" ? NULL : periodicFile.c_str(),
					      check);

      // If there is no value file which should be written, we can delete
      // the information of the macro file.
      if (!valueFilename.length()) {
	clearMacroFileInfo();
      }
    }

    initialized = true;
  }

  bool Mesh::associated(DegreeOfFreedom dof1, DegreeOfFreedom dof2) {
    ::std::map<BoundaryType, VertexVector*>::iterator it;
    ::std::map<BoundaryType, VertexVector*>::iterator end = periodicAssociations.end();
    for (it = periodicAssociations.begin(); it != end; ++it) {
      if ((*(it->second))[dof1] == dof2)
	return true;
    }
    return false;
  }

  bool Mesh::indirectlyAssociated(DegreeOfFreedom dof1, DegreeOfFreedom dof2) {
    ::std::vector<DegreeOfFreedom> associatedToDOF1;
    int i, size;
    ::std::map<BoundaryType, VertexVector*>::iterator it;
    ::std::map<BoundaryType, VertexVector*>::iterator end = periodicAssociations.end();
    DegreeOfFreedom dof, assDOF;

    associatedToDOF1.push_back(dof1);
    for(it = periodicAssociations.begin(); it != end; ++it) {
      size = static_cast<int>(associatedToDOF1.size());
      for(i = 0; i < size; i++) {
	dof = associatedToDOF1[i];
	assDOF = (*(it->second))[dof];
	if(assDOF == dof2) {
	  return true;
	} else {
	  if(assDOF != dof) associatedToDOF1.push_back(assDOF);
	}
      }
    }
    return false;
  }

  void Mesh::clearMacroFileInfo()
  {
    macroFileInfo_->clear(getNumberOfEdges(),
			  getNumberOfVertices());
    DELETE macroFileInfo_;
  }
}