void BaseMoveableMem :: MoveTo (size_t newpos)
  {
    //  cout << "move block, oldpos = " << pos << "; newpos = " << newpos
    // << ", size = " << size << endl;
    static size_t move = 0;

    if (newpos + size > totalsize)
      throw NgException ("MoveableMem overflow");
    if (newpos > pos)
      {
        if (next) next->MoveTo (newpos+size);
        memmove (largeblock+newpos, largeblock+pos, size);
        move += size;
      }
    else if (newpos < pos)
      {
        // cout << "move down: " << size << endl;
        memmove (largeblock+newpos, largeblock+pos, size);
        if (next) next->MoveTo (newpos+size);
        move += size;
      }
    pos = newpos;
    ptr = largeblock+pos;
    //  cout << "total move: " << move << endl;
  }
void SplineGeometry2d :: Load (const char * filename)
{
    ifstream infile;
    Point<2> x;
    char buf[50];

    infile.open (filename);

    if ( ! infile.good() )
        throw NgException(string ("Input file '") +
                          string (filename) +
                          string ("' not available!"));

    TestComment ( infile );

    infile >> buf;   // file recognition

    tensormeshing.SetSize(0);
    quadmeshing.SetSize(0);

    TestComment ( infile );
    if ( strcmp (buf, "splinecurves2dnew") == 0 )
    {
        LoadDataNew ( infile );
    }
    else if ( strcmp (buf, "splinecurves2dv2") == 0 )
    {
        LoadDataV2 ( infile );
    }
    else
    {
        LoadData(infile );
    }
    infile.close();
}
Beispiel #3
0
  bool BASE_INDEX_3_CLOSED_HASHTABLE ::
  PositionCreate2 (const INDEX_3 & ind, int & apos) 
  {
    int i = HashValue(ind);
    int startpos = i;
    while (1)
      {
        /*
	i++;
	if (i >= hash.Size()) i = 0;
        */
        i = (i+1) % hash.Size();
	if (hash[i] == ind) 
	  {
	    apos = i;
	    return false;
	  }
	if (hash[i].I1() == invalid) 
	  {
	    hash[i] = ind;
	    apos = i;
	    return true;
	  }
	if (i == startpos)
	  throw NgException ("Try to set new element in full closed hashtable");
      }
  }
Beispiel #4
0
  void GeomSearch3d :: AddElem(const MiniElement2d& elem, INDEX elemnum)
  {
    Point3d minp, maxp;
    ElemMaxExt(minp, maxp, elem);
    int sx = int ((minp.X()-minext.X())/elemsize.X()+1.);
    int ex = int ((maxp.X()-minext.X())/elemsize.X()+1.);
    int sy = int ((minp.Y()-minext.Y())/elemsize.Y()+1.);
    int ey = int ((maxp.Y()-minext.Y())/elemsize.Y()+1.);
    int sz = int ((minp.Z()-minext.Z())/elemsize.Z()+1.);
    int ez = int ((maxp.Z()-minext.Z())/elemsize.Z()+1.);
  
    for (int ix = sx; ix <= ex; ix++)
      for (int iy = sy; iy <= ey; iy++)
        for (int iz = sz; iz <= ez; iz++)
          {
            INDEX ind=ix+(iy-1)*size.i1+(iz-1)*size.i2*size.i1;
            if (ind < 1 || ind > size.i1 * size.i2 * size.i3)
              {
                cerr << "Illegal hash-position";
                cerr << "Position: " << ix << "," << iy << "," << iz << endl;
		    throw NgException ("Illegal position in Geomsearch");
              }
            hashtable.Elem(ind)->Append(elemnum);		      
          }
  }
Beispiel #5
0
int Polyhedra :: AddFace (int pi1, int pi2, int pi3, int inputnum)
{
  (*testout) << "polyhedra, add face " << pi1 << ", " << pi2 << ", " << pi3 << endl;

  if(pi1 == pi2 || pi2 == pi3 || pi3 == pi1)
    {
      ostringstream msg;
      msg << "Illegal point numbers for polyhedron face: " << pi1+1 << ", " << pi2+1 << ", " << pi3+1;
      throw NgException(msg.str());
    }

  faces.Append (Face (pi1, pi2, pi3, points, inputnum));
  
  Point<3> p1 = points[pi1];
  Point<3> p2 = points[pi2];
  Point<3> p3 = points[pi3];

  Vec<3> v1 = p2 - p1;
  Vec<3> v2 = p3 - p1;

  Vec<3> n = Cross (v1, v2); 
  n.Normalize();

  Plane pl (p1, n);
//   int inverse;
//   int identicto = -1;
//   for (int i = 0; i < planes.Size(); i++)
//     if (pl.IsIdentic (*planes[i], inverse, 1e-9*max3(v1.Length(),v2.Length(),Dist(p2,p3))))
//       {
// 	if (!inverse)
// 	  identicto = i;
//       }
//   //  cout << "is identic = " << identicto << endl;
//   identicto = -1;    // changed April 10, JS

//   if (identicto != -1)
//     faces.Last().planenr = identicto;
//   else
    {
      planes.Append (new Plane (p1, n));
      surfaceactive.Append (1);
      surfaceids.Append (0);
      faces.Last().planenr = planes.Size()-1;
    }

//  (*testout) << "is plane nr " << faces.Last().planenr << endl;

  return faces.Size();
}
Beispiel #6
0
Primitive * Primitive :: CreatePrimitive (const char * classname)
{
  if (strcmp (classname, "sphere") == 0)
    return Sphere::CreateDefault();
  if (strcmp (classname, "plane") == 0)
    return Plane::CreateDefault();
  if (strcmp (classname, "cylinder") == 0)
    return Cylinder::CreateDefault();
  if (strcmp (classname, "cone") == 0)
    return Cone::CreateDefault();
  if (strcmp (classname, "brick") == 0)
    return Brick::CreateDefault();


  stringstream ost;
  ost << "Primitve::CreatePrimitive not implemented for " << classname << endl;
  throw NgException (ost.str());
}
  int BASE_INDEX_CLOSED_HASHTABLE ::
  PositionCreate2 (const INDEX & ind, int & apos) 
  {
    int i = HashValue(ind);
    int startpos = i;
    while (1)
      {
	i++;
	if (i > hash.Size()) i = 1;
	if (hash.Get(i) == ind) 
	  {
	    apos = i;
	    return 0;
	  }
	if (hash.Get(i) == invalid) 
	  {
	    hash.Elem(i) = ind;
	    apos = i;
	    return 1;
	  }
	if (i == startpos)
	  throw NgException ("Try to set new element in full closed hashtable");
      }
  }
Beispiel #8
0
  void WriteTETFormat (const Mesh & mesh,
		       const string & filename)//, const string& problemType )
  {
    string problemType = "";
    if(!mesh.PureTetMesh())
      throw NgException("Can only export pure tet mesh in this format");

    cout << "starting .tet export to file " << filename << endl;


    ARRAY<int> point_ids,edge_ids,face_ids;
    ARRAY<int> elnum(mesh.GetNE());
    elnum = -1;

    
    ARRAY<int> userdata_int;
    ARRAY<double> userdata_double;
    ARRAY<int> ports;

    ARRAY<int> uid_to_group_3D, uid_to_group_2D, uid_to_group_1D, uid_to_group_0D;

    int pos_int = 0;
    int pos_double = 0;
    
    bool haveuserdata = 
      (mesh.GetUserData("TETmesh:double",userdata_double) &&
       mesh.GetUserData("TETmesh:int",userdata_int) && 
       mesh.GetUserData("TETmesh:ports",ports) &&
       mesh.GetUserData("TETmesh:point_id",point_ids,PointIndex::BASE) &&
       mesh.GetUserData("TETmesh:uid_to_group_3D",uid_to_group_3D) &&
       mesh.GetUserData("TETmesh:uid_to_group_2D",uid_to_group_2D) &&
       mesh.GetUserData("TETmesh:uid_to_group_1D",uid_to_group_1D) &&
       mesh.GetUserData("TETmesh:uid_to_group_0D",uid_to_group_0D));


    int version,subversion;

    if(haveuserdata)
      {
	version = int(userdata_double[0]);
	subversion = int(10*(userdata_double[0] - version));
	pos_double++;
      }
    else
      {
	version = 2;
	subversion = 0;
      }

    
    if(version >= 2)
      {
	// test if ids are disjunct, if not version 2.0 not possible
	int maxbc(-1),mindomain(-1);
	
	for(ElementIndex i=0; i<mesh.GetNE(); i++)
	  if(i==0 || mesh[i].GetIndex() < mindomain)
	    mindomain = mesh[i].GetIndex();
	for(int i=1; i<=mesh.GetNFD(); i++)
	  if(i==1 || mesh.GetFaceDescriptor(i).BCProperty() > maxbc)
	    maxbc = mesh.GetFaceDescriptor(i).BCProperty();
	
	if(maxbc >= mindomain)
	  {
	    cout << "WARNING: writing version " << version << "." << subversion << " tetfile not possible, ";
	    version = 1; subversion = 1;
	    cout << "using version " << version << "." << subversion << endl;
	  }
      }



    int startsize = point_ids.Size();
    point_ids.SetSize(mesh.GetNP()+1);
    for(int i=startsize; i<point_ids.Size(); i++)
      point_ids[i] = -1;


    for(int i=0; i<PointIndex::BASE; i++)
      point_ids[i] = -1;


    INDEX_2_CLOSED_HASHTABLE<int> edgenumbers(6*mesh.GetNE()+3*mesh.GetNSE());;
    INDEX_3_CLOSED_HASHTABLE<int> facenumbers(4*mesh.GetNE()+mesh.GetNSE());

    ARRAY<INDEX_2> edge2node;
    ARRAY<INDEX_3> face2edge;
    ARRAY<INDEX_4> element2face;

    int numelems(0),numfaces(0),numedges(0),numnodes(0);

    for(SegmentIndex si = 0; si < mesh.GetNSeg(); si++)
      {
	const Segment & seg = mesh[si];
	INDEX_2 i2(seg.p1,seg.p2);
	i2.Sort();
	if(edgenumbers.Used(i2))
	  continue;

	numedges++;
	edgenumbers.Set(i2,numedges);
	edge2node.Append(i2);

	edge_ids.Append(seg.edgenr);

	if(point_ids[seg.p1] == -1)
	  point_ids[seg.p1] = (version >= 2) ? seg.edgenr : 0;
	if(point_ids[seg.p2] == -1)
	  point_ids[seg.p2] = (version >= 2) ? seg.edgenr : 0;
      }

    for(SurfaceElementIndex si = 0; si < mesh.GetNSE(); si++)
      {
	if(mesh[si].IsDeleted())
	  continue;

	const Element2d & elem = mesh[si];

	numfaces++;
	INDEX_3 i3(elem[0], elem[1], elem[2]);

	int min = i3[0];
	int minpos = 0;
	for(int j=1; j<3; j++)
	  if(i3[j] < min)
	    {
	      min = i3[j]; minpos = j;
	    }
	if(minpos == 1)
	  {
	    int aux = i3[0]; i3[0] = i3[1]; i3[1] = i3[2]; i3[2] = aux;
	  }
	else if(minpos == 2)
	  {
	    int aux = i3[0]; i3[0] = i3[2]; i3[2] = i3[1]; i3[1] = aux;
	  }
	facenumbers.Set(i3,numfaces);

	int bc = mesh.GetFaceDescriptor(elem.GetIndex()).BCProperty();
	face_ids.Append(bc);

	for(int j=0; j<3; j++)
	  if(point_ids[elem[j]] == -1)
	    point_ids[elem[j]] = (version >= 2) ? bc : 0;

	INDEX_2 i2a,i2b;
	INDEX_3 f_to_n;
	for(int j=0; j<3; j++)
	  {
	    i2a = INDEX_2(i3[j],i3[(j+1)%3]);
	    i2b[0] = i2a[1]; i2b[1] = i2a[0];
	    if(edgenumbers.Used(i2a))
	      f_to_n[j] = edgenumbers.Get(i2a);
	    else if(edgenumbers.Used(i2b))
	      f_to_n[j] = -edgenumbers.Get(i2b);
	    else
	      {
		numedges++;
		edgenumbers.Set(i2a,numedges);
		edge2node.Append(i2a);
		f_to_n[j] = numedges;
		if(version >= 2)
		  edge_ids.Append(bc);
		else
		  edge_ids.Append(0);
	      }
	  }
	face2edge.Append(f_to_n);
      }
    
    for(ElementIndex ei = 0; ei < mesh.GetNE(); ei++)
      {
	const Element & el = mesh[ei];

	if(el.IsDeleted())
	  continue;

	numelems++;
	elnum[ei] = numelems;

	static int tetfaces[4][3] =
	  { { 0, 2, 1 },
	    { 0, 1, 3 },
	    { 1, 2, 3 },
	    { 2, 0, 3 } };
	
	for(int j=0; j<4; j++)
	  if(point_ids[el[j]] == -1)
	    point_ids[el[j]] = (version >= 2) ? el.GetIndex() : 0;

	INDEX_4 e_to_f;

	for(int i = 0; i < 4; i++)
	  {
	    INDEX_3 i3a(el[tetfaces[i][0]],el[tetfaces[i][1]],el[tetfaces[i][2]]);
	    
	    int min = i3a[0];
	    int minpos = 0;
	    for(int j=1; j<3; j++)
	      if(i3a[j] < min)
		{
		  min = i3a[j]; minpos = j;
		}
	    if(minpos == 1)
	      {
		int aux = i3a[0]; i3a[0] = i3a[1]; i3a[1] = i3a[2]; i3a[2] = aux;
	      }
	    else if(minpos == 2)
	      {
		int aux = i3a[0]; i3a[0] = i3a[2]; i3a[2] = i3a[1]; i3a[1] = aux;
	      }
	    INDEX_3 i3b(i3a[0],i3a[2],i3a[1]);
	    

	    if(facenumbers.Used(i3a))
	      e_to_f[i] = facenumbers.Get(i3a);
	    else if(facenumbers.Used(i3b))
	      e_to_f[i] = -facenumbers.Get(i3b);
	    else
	      {
		numfaces++;
		facenumbers.Set(i3a,numfaces);
		e_to_f[i] = numfaces;
		if(version >= 2)
		  face_ids.Append(el.GetIndex());
		else
		  face_ids.Append(0);

		INDEX_2 i2a,i2b;
		INDEX_3 f_to_n;
		for(int j=0; j<3; j++)
		  {
		    i2a = INDEX_2(i3a[j],i3a[(j+1)%3]);
		    i2b[0] = i2a[1]; i2b[1] = i2a[0];
		    if(edgenumbers.Used(i2a))
		      f_to_n[j] = edgenumbers.Get(i2a);
		    else if(edgenumbers.Used(i2b))
		      f_to_n[j] = -edgenumbers.Get(i2b);
		    else
		      {
			numedges++;
			edgenumbers.Set(i2a,numedges);
			edge2node.Append(i2a);
			f_to_n[j] = numedges;
			if(version >= 2)
			  edge_ids.Append(el.GetIndex());
			else
			  edge_ids.Append(0);
		      }
		  }
		face2edge.Append(f_to_n);	  
	      }
	  }
	element2face.Append(e_to_f);
      }




    ofstream outfile(filename.c_str());

    outfile.precision(16);

    int unitcode;
    double tolerance;
    double dS1,dS2, alphaDeg;
    double x3D,y3D,z3D;
    int modelverts(0), modeledges(0), modelfaces(0), modelcells(0);

    int numObj0D,numObj1D,numObj2D,numObj3D;
    int numports = ports.Size();

    ARRAY<int> nodenum(point_ids.Size()+1);

    nodenum = -1;
	    


    numnodes = 0;
    for(int i=0; i<point_ids.Size(); i++)
      {
	if(point_ids[i] != -1)
	  {
	    numnodes++;
	    nodenum[i] = numnodes;
	  }
      }


    if(haveuserdata)
      {
	unitcode = userdata_int[pos_int];
	pos_int++;

	tolerance = userdata_double[pos_double];
	pos_double++;

	dS1 = userdata_double[pos_double];
	pos_double++;
	dS2 = userdata_double[pos_double];
	pos_double++;
	alphaDeg = userdata_double[pos_double];
	pos_double++;

	x3D = userdata_double[pos_double];
	pos_double++;
	y3D = userdata_double[pos_double];
	pos_double++;
	z3D = userdata_double[pos_double];
	pos_double++;

	if(version == 2)
	  {
	    modelverts = userdata_int[pos_int];
	    pos_int++;
	    modeledges = userdata_int[pos_int];
	    pos_int++;
	    modelfaces = userdata_int[pos_int];
	    pos_int++;
	    modelcells = userdata_int[pos_int];
	    pos_int++;
	  }

	numObj3D = userdata_int[pos_int];
	pos_int++;
	numObj2D = userdata_int[pos_int];
	pos_int++;
	numObj1D = userdata_int[pos_int];
	pos_int++;
	numObj0D = userdata_int[pos_int];
	pos_int++;
      }
    else
      {
	unitcode = 3;

	tolerance = 1e-5;

	dS1 = dS2 = alphaDeg = 0;

	x3D = y3D = z3D = 0;

	modelverts = modeledges = modelfaces = modelcells = 0;
	
	numObj3D = numObj2D = numObj1D = numObj0D = 0;
      }

    string uidpid;
    if(version == 1)
      uidpid = "PID";
    else if (version == 2)
      uidpid = "UID";
    

    ARRAY< ARRAY<int,PointIndex::BASE>* > idmaps;
    for(int i=1; i<=mesh.GetIdentifications().GetMaxNr(); i++)
      {
	if(mesh.GetIdentifications().GetType(i) == Identifications::PERIODIC)
	  {
	    idmaps.Append(new ARRAY<int,PointIndex::BASE>);
	    mesh.GetIdentifications().GetMap(i,*idmaps.Last(),true);
	  }
      }

    ARRAY<int> id_num,id_type;
    ARRAY< ARRAY<int> *> id_groups;


	// sst 2008-03-12: Write problem class...
	{
		std::string block;
		block  = "// CST Tetrahedral ";
		block += !problemType.empty() ? problemType : "High Frequency";
		block += " Mesh, Version no.:\n";
		
		size_t size = block.size()-3;
		block += "// ";
		block.append( size, '^' );
		block += "\n";

		outfile
			<< block
			<< version << "." << subversion << "\n\n";
	}

	outfile 
	    << "// User Units Code (1=CM 2=MM 3=M 4=MIC 5=NM 6=FT 7=IN 8=MIL):\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n" \
	    << unitcode << "\n\n"					\
	    << "// Geometric coord \"zero\" tolerance threshold:\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n" \
	    << tolerance << "\n\n"				  \
	    << "// Periodic UnitCell dS1 , dS2 , alphaDeg:\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n" \
	    << dS1 << " " << dS2 << " " << alphaDeg <<"\n\n"	\
	    << "// Periodic UnitCell origin in global coords (x3D,y3D,z3D):\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n" \
	    << x3D << " " << y3D << " " << z3D << "\n" << endl;

    if(version == 2)
      {
	outfile << "// Model entity count: Vertices, Edges, Faces, Cells:\n" \
		<< "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n" \
		<< modelverts << " " << modeledges << " " << modelfaces << " " << modelcells << endl << endl;
      }


    outfile << "// Topological mesh-entity counts (#elements,#faces,#edges,#nodes):\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n";
    outfile << numelems << " " 
	    << numfaces << " "
	    << numedges << " " 
	    << numnodes << endl << endl;

    outfile << "// NodeID, X, Y, Z, Type (0=Reg 1=PMaster 2=PSlave 3=CPMaster 4=CPSlave), "<< uidpid <<":\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n";
       


    id_num.SetSize(mesh.GetNP()+1);
    id_type.SetSize(mesh.GetNP()+1);
    id_num = 0;
    id_type = 0;

    int n2,n4,n8;
    n2 = n4 = n8 = 0;

 
    for(int i=PointIndex::BASE; i<mesh.GetNP()+PointIndex::BASE; i++)
      {
	if(id_num[i] != 0)
	  continue;

	if(nodenum[i] == -1)
	  continue;

	ARRAY<int> group;
	group.Append(i);
	for(int j=0; j<idmaps.Size(); j++)
	  {
	    startsize = group.Size();
	    for(int k=0; k<startsize; k++)
	      {
		int id = (*idmaps[j])[group[k]];
		if(id != 0 && !group.Contains(id) && nodenum[id] != -1)
		  {
		    group.Append(id);
		    id_num[id] = j+1+id_num[group[k]];
		  }
	      }
	  }
	if(group.Size() > 1)
	  {
	    id_groups.Append(new ARRAY<int>(group));
	    if(group.Size() == 2)
	      {
		id_type[i] = 1;
		id_type[group[1]] = 2;
		n2++;
	      }
	    else if(group.Size() == 4)
	      {
		id_type[i] = 3;
		for(int j=1; j<group.Size(); j++)
		  id_type[group[j]] = 4;
		n4++;
	      }
	    else if(group.Size() == 8)
	      {
		id_type[i] = 5;
		for(int j=1; j<group.Size(); j++)
		  id_type[group[j]] = 6;
		n8++;
	      }
	    else
	      cerr << "ERROR: Identification group size = " << group.Size() << endl;
	  }
	
      }


    for(PointIndex i=PointIndex::BASE; i<mesh.GetNP()+PointIndex::BASE; i++)
      {
	if(nodenum[i] == -1)
	  continue;
	outfile << nodenum[i] << " "
		<< mesh[i](0) << " "
		<< mesh[i](1) << " "
		<< mesh[i](2) << " " << id_type[i] << " ";
	if(i-PointIndex::BASE < point_ids.Size())
	  outfile << point_ids[i];
	else
	  outfile << "0";
	outfile << "\n";
      }
    outfile << endl;

    outfile << "\n// Number of Periodic Master Nodes:\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n" \
	    << n2 << "\n"			       \
	    << "\n" \
	    << "// MasterNodeID, SlaveNodeID, TranslCode (1=dS1 2=dS2 3=dS1+dS2):\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n";
    for(int i=0; i<id_groups.Size(); i++)
      {
	if(id_groups[i]->Size() != 2)
	  continue;

	for(int j=0; j<id_groups[i]->Size(); j++)
	  outfile << nodenum[(*id_groups[i])[j]] << " ";
	for(int j=1; j<id_groups[i]->Size(); j++)
	  outfile << id_num[(*id_groups[i])[j]] << " ";
	outfile << "\n";

	delete id_groups[i];
	id_groups[i] = NULL;
      }
    outfile << endl;
	
	
    outfile << "// Number of Corner Periodic Master Nodes:\n"	      \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n" \
	    << n4 << "\n"				      \
	    << "\n" \
	    << "// MasterNodeID, 3-SlaveNodeID's, 3-TranslCodes (1=dS1 2=dS2 3=dS1+dS2):\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n";


    for(int i=0; i<id_groups.Size(); i++)
      {
	if(!id_groups[i] || id_groups[i]->Size() != 4)
	  continue;

	for(int j=0; j<id_groups[i]->Size(); j++)
	  outfile << nodenum[(*id_groups[i])[j]] << " ";
	for(int j=1; j<id_groups[i]->Size(); j++)
	  {
	    outfile << id_num[(*id_groups[i])[j]] << " ";
	  }
	outfile << "\n";

	delete id_groups[i];
	id_groups[i] = NULL;
      }
    outfile << endl;


    outfile << "// Number of Cubic Periodic Master Nodes:\n"	     \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n" \
	    << n8 << "\n"				     \
	    << "\n" \
	    << "// MasterNodeID, 7-SlaveNodeID's, TranslCodes:\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n";
    for(int i=0; i<id_groups.Size(); i++)
      {
	if(!id_groups[i] || id_groups[i]->Size() != 8)
	  continue;

	for(int j=0; j<id_groups[i]->Size(); j++)
	  outfile << nodenum[(*id_groups[i])[j]] << " ";
	for(int j=1; j<id_groups[i]->Size(); j++)
	  outfile << id_num[(*id_groups[i])[j]] << " ";
	outfile << "\n";

	delete id_groups[i];
	id_groups[i] = NULL;
      }
    outfile << endl;

    


    outfile << "// EdgeID, NodeID0, NodeID1, Type (0=Reg 1=PMaster 2=PSlave 3=CPMaster 4=CPSlave), "<<uidpid<<":\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n";

    
      
    ARRAY< ARRAY<int>* > vertex_to_edge(mesh.GetNP()+1);
    for(int i=0; i<=mesh.GetNP(); i++)
      vertex_to_edge[i] = new ARRAY<int>;

    ARRAY< ARRAY<int,PointIndex::BASE>* > idmaps_edge(idmaps.Size());
    for(int i=0; i<idmaps_edge.Size(); i++)
      {
	idmaps_edge[i] = new ARRAY<int,PointIndex::BASE>(numedges);
	(*idmaps_edge[i]) = 0;
      }

    ARRAY<int> possible;
    for(int i=0; i<edge2node.Size(); i++)
      {
	const INDEX_2 & v = edge2node[i];
	for(int j=0; j<idmaps.Size(); j++)
	  {
	    INDEX_2 vid((*idmaps[j])[v[0]], (*idmaps[j])[v[1]]);
	    if(vid[0] != 0 && vid[0] != v[0] && vid[1] != 0 && vid[1] != v[1])
	      {
		Intersection(*vertex_to_edge[vid[0]],*vertex_to_edge[vid[1]],possible);
		if(possible.Size() == 1)
		  {
		    (*idmaps_edge[j])[possible[0]] = i+1;
		    (*idmaps_edge[j])[i+1] = possible[0];
		  }
		else if(possible.Size() > 0)
		  {
		    cerr << "ERROR: too many possible edge identifications" << endl;
		    (*testout) << "ERROR: too many possible edge identifications" << endl
			       << "*vertex_to_edge["<<vid[0]<<"] " << *vertex_to_edge[vid[0]] << endl
			       << "*vertex_to_edge["<<vid[1]<<"] " << *vertex_to_edge[vid[1]] << endl
			       << "possible " << possible << endl;
		  }
	      }
	  }
	vertex_to_edge[v[0]]->Append(i+1);
	vertex_to_edge[v[1]]->Append(i+1);
      }


    for(int i=0; i<vertex_to_edge.Size(); i++)
      delete vertex_to_edge[i];


    id_groups.SetSize(0);
    id_num.SetSize(numedges+1);
    id_num = 0;
    id_type.SetSize(numedges+1);
    id_type = 0;

    n2 = n4 = n8 = 0;

    for(int i=1; i<=edge2node.Size(); i++)
      {
	if(id_num[i] != 0)
	  continue;


	ARRAY<int> group;
	group.Append(i);
	for(int j=0; j<idmaps_edge.Size(); j++)
	  {
	    startsize = group.Size();
	    for(int k=0; k<startsize; k++)
	      {
		int id = (*idmaps_edge[j])[group[k]];
		if(id != 0 && !group.Contains(id))
		  {
		    group.Append(id);
		    id_num[id] = j+1+id_num[group[k]];
		  }
	      }
	  }
	if(group.Size() > 1)
	  {
	    id_num[i] = 1;
	    id_groups.Append(new ARRAY<int>(group));
	    if(group.Size() == 2)
	      {
		id_type[i] = 1;
		id_type[group[1]] = 2;
		n2++;
	      }
	    else if(group.Size() == 4)
	      {
		id_type[i] = 3;
		for(int j=1; j<group.Size(); j++)
		  id_type[group[j]] = 4;
		n4++;
	      }
	    else
	      {
		cerr << "ERROR: edge identification group size = " << group.Size() << endl;
		(*testout) << "edge group " << group << endl;
		for(int j=0; j<idmaps_edge.Size(); j++)
		  {
		    (*testout) << "edge id map " << j << endl << *idmaps_edge[j] << endl;
		  }
	      }
	  }
      }



    for(int i=1; i<=edge2node.Size(); i++)
      {
	if(id_num[i] != 0)
	  continue;


	ARRAY<int> group;
	group.Append(i);
	for(int j=0; j<idmaps_edge.Size(); j++)
	  {
	    startsize = group.Size();
	    for(int k=0; k<startsize; k++)
	      {
		int id = (*idmaps_edge[j])[group[k]];
		if(id != 0 && !group.Contains(id))
		  {
		    group.Append(id);
		    id_num[id] = j+1+id_num[group[k]];
		  }
	      }
	  }
	if(group.Size() > 1)
	  {
	    id_num[i] = 1;
	    id_groups.Append(new ARRAY<int>(group));
	    if(group.Size() == 2)
	      {
		id_type[i] = 1;
		id_type[group[1]] = 2;
		n2++;
	      }
	    else if(group.Size() == 4)
	      {
		id_type[i] = 3;
		for(int j=1; j<group.Size(); j++)
		  id_type[group[j]] = 4;
		n4++;
	      }
	    else
	      {
		cerr << "ERROR: edge identification group size = " << group.Size() << endl;
		(*testout) << "edge group " << group << endl;
		for(int j=0; j<idmaps_edge.Size(); j++)
		  {
		    (*testout) << "edge id map " << j << endl << *idmaps_edge[j] << endl;
		  }
	      }
	  }
	
      }

    
    for(int i=0; i<edge2node.Size(); i++)
      outfile << i+1 << " " << nodenum[edge2node[i][0]] << " " << nodenum[edge2node[i][1]] 
	      << " " << id_type[i+1] << " " << edge_ids[i] << "\n";

    outfile << endl;

    

    outfile << "// Number of Periodic Master Edges:\n"\
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n"\
	    << n2 << "\n"			      \
	    << "\n"\
	    << "// MasterEdgeID, SlaveEdgeID, TranslCode (1=dS1 2=dS2 3=dS1+dS2):\n"\
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n";
    for(int i=0; i<id_groups.Size(); i++)
      {
	if(id_groups[i]->Size() != 2)
	  continue;

	for(int j=0; j<id_groups[i]->Size(); j++)
	  outfile << (*id_groups[i])[j] << " ";
	for(int j=1; j<id_groups[i]->Size(); j++)
	  outfile << id_num[(*id_groups[i])[j]] << " ";
	outfile << "\n";

	delete id_groups[i];
	id_groups[i] = NULL;
      }
    outfile << endl;

    outfile << "// Number of Corner Periodic Master Edges:\n"		\
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n"\
	    << n4 << "\n"				     \
	    << "\n"\
	    << "// MasterEdgeID, 3 SlaveEdgeID's, 3 TranslCode (1=dS1 2=dS2 3=dS1+dS2):\n"\
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n";
    for(int i=0; i<id_groups.Size(); i++)
      {
	if(!id_groups[i] || id_groups[i]->Size() != 4)
	  continue;

	for(int j=0; j<id_groups[i]->Size(); j++)
	  outfile << (*id_groups[i])[j] << " ";
	for(int j=1; j<id_groups[i]->Size(); j++)
	  outfile << id_num[(*id_groups[i])[j]] << " ";
	outfile << "\n";

	delete id_groups[i];
	id_groups[i] = NULL;
      }
    outfile << endl;


    outfile << "// FaceID, EdgeID0, EdgeID1, EdgeID2, FaceType (0=Reg 1=PMaster 2=PSlave), "<<uidpid<<":\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n";

    
    
    ARRAY< ARRAY<int>* > edge_to_face(numedges+1);
    for(int i=0; i<edge_to_face.Size(); i++)
      edge_to_face[i] = new ARRAY<int>;

    
    for(int i=0; i<idmaps.Size(); i++)
      {
	idmaps[i]->SetSize(numfaces);
	(*idmaps[i]) = 0;
      }

    
    for(int i=0; i<face2edge.Size(); i++)
      {
	for(int j=0; j<idmaps_edge.Size(); j++)
	  {
	    int e1id,e2id,e3id;
	    e1id = (*idmaps_edge[j])[abs(face2edge[i][0])];
	    e2id = (*idmaps_edge[j])[abs(face2edge[i][1])];
	    e3id = (*idmaps_edge[j])[abs(face2edge[i][2])];
	    if(e1id != 0 && e1id != abs(face2edge[i][0]) &&
	       e2id != 0 && e2id != abs(face2edge[i][1]) &&
	       e3id != 0 && e3id != abs(face2edge[i][2]))
	      {
		Intersection(*edge_to_face[e1id],*edge_to_face[e2id],*edge_to_face[e3id],possible);
		if(possible.Size() == 1)
		  {
		    (*idmaps[j])[possible[0]] = i+1;
		    (*idmaps[j])[i+1] = possible[0];
		  }
		else if(possible.Size() > 0)
		  cerr << "ERROR: too many possible face identifications" << endl;
	      }
	  }

	edge_to_face[abs(face2edge[i][0])]->Append(i+1);
	edge_to_face[abs(face2edge[i][1])]->Append(i+1);
	edge_to_face[abs(face2edge[i][2])]->Append(i+1);
      }

    for(int i=0; i<edge_to_face.Size(); i++)
      delete edge_to_face[i];


    for(int i=0; i<idmaps_edge.Size(); i++)
      delete idmaps_edge[i];

    
    id_groups.SetSize(0);
    id_num.SetSize(numfaces+1);
    id_num = 0;

    n2 = n4 = n8 = 0;

    for(int i=1; i<=numfaces; i++)
      {
	if(id_num[i] != 0)
	  continue;

	ARRAY<int> group;
	group.Append(i);
	for(int j=0; j<idmaps.Size(); j++)
	  {
	    startsize = group.Size();
	    for(int k=0; k<startsize; k++)
	      {
		int id = (*idmaps[j])[group[k]];
		if(id != 0 && !group.Contains(id))
		  {
		    group.Append(id);
		    id_num[id] = j+1+id_num[group[k]];
		  }
	      }
	  }
	if(group.Size() > 1)
	  {
	    id_num[i] = -1;
	    id_groups.Append(new ARRAY<int>(group));
	    if(group.Size() == 2)
	      n2++;
	    else
	      cerr << "ERROR: face identification group size = " << group.Size() << endl;
	  }
	
      }


    for(int i=0; i<idmaps.Size(); i++)
      delete idmaps[i];




    for(int i=0; i<face2edge.Size(); i++)
      {	
	outfile << i+1 << " ";
	for(int j=0; j<3; j++)
	  outfile << face2edge[i][j] << " ";

	if(id_num[i+1] == 0)
	  outfile << 0;
	else if(id_num[i+1] == -1)
	  outfile << 1;
	else
	  outfile << 2;

	outfile << " " << face_ids[i] <<"\n";
      }
    outfile << endl;


    outfile << "// Number of Periodic Master Faces:\n"\
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n"\
	    << n2 << "\n"			      \
	    << "\n"\
	    << "// MasterFaceID, SlaveFaceID, TranslCode (1=dS1 2=dS2):\n"\
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n";
    for(int i=0; i<id_groups.Size(); i++)
      {
	if(id_groups[i]->Size() != 2)
	  continue;

	for(int j=0; j<id_groups[i]->Size(); j++)
	  outfile << (*id_groups[i])[j] << " ";
	for(int j=1; j<id_groups[i]->Size(); j++)
	  outfile << id_num[(*id_groups[i])[j]] << " ";
	outfile << "\n";

	delete id_groups[i];
      }
    outfile << endl;

    


    outfile << "// ElemID, FaceID0, FaceID1, FaceID2, FaceID3, "<<uidpid<<":\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n";

    for(ElementIndex i=0; i<mesh.GetNE(); i++)
      {
	if(elnum[i] >= 0)
	  {
	    outfile << elnum[i] << " ";
	    for(int j=0; j<4; j++)
	      outfile << element2face[elnum[i]-1][j] << " ";

	    outfile << mesh[i].GetIndex() << "\n";
	  }
      }
    outfile << endl;

    outfile << "// ElemID, NodeID0, NodeID1, NodeID2, NodeID3:\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n";

    
    for(ElementIndex i=0; i<mesh.GetNE(); i++)
      {
	if(elnum[i] >= 0)
	  outfile << elnum[i] << " "
		  << nodenum[mesh[i][1]] << " " << nodenum[mesh[i][0]] << " " << nodenum[mesh[i][2]] << " " << nodenum[mesh[i][3]] << "\n";
      }
    outfile << endl;
    

    

    outfile << "// Physical Object counts (#Obj3D,#Obj2D,#Obj1D,#Obj0D):\n"
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n"
	    << " "<< numObj3D << " " << numObj2D << " " << numObj1D << " " << numObj0D << "\n" \
	    << "\n" \
	    << "// Number of Ports (Ports are a subset of Object2D list):\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n" \
	    << numports << "\n"						\
	    << endl;


    ARRAY< ARRAY<int> * > groups;

    int maxg = -1;
    for(int i = 0; i<uid_to_group_3D.Size(); i++)
      if(uid_to_group_3D[i] > maxg)
	maxg = uid_to_group_3D[i];
    for(int i = 0; i<uid_to_group_2D.Size(); i++)
      if(uid_to_group_2D[i] > maxg)
	maxg = uid_to_group_2D[i];
    for(int i = 0; i<uid_to_group_1D.Size(); i++)
      if(uid_to_group_1D[i] > maxg)
	maxg = uid_to_group_1D[i];
    for(int i = 0; i<uid_to_group_0D.Size(); i++)
      if(uid_to_group_0D[i] > maxg)
	maxg = uid_to_group_0D[i];

    groups.SetSize(maxg+1);
    for(int i=0; i<groups.Size(); i++)
      groups[i] = new ARRAY<int>;

    for(ElementIndex i=0; i<mesh.GetNE(); i++)
      if(uid_to_group_3D[mesh[i].GetIndex()] >= 0)
	groups[uid_to_group_3D[mesh[i].GetIndex()]]->Append(i+1);
      
    


    outfile << "// Object3D GroupID, #Elems <immediately followed by> ElemID List:\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n";
    for(int i=0; i<numObj3D; i++)
      {
	outfile << i << " " << groups[i]->Size() << "\n";
	for(int j=0; j<groups[i]->Size(); j++)
	  outfile << (*groups[i])[j] << "\n";
      }

    for(int i=0; i<groups.Size(); i++)
      groups[i]->SetSize(0);

    for(int i=0; i<face_ids.Size(); i++)
      if(uid_to_group_2D[face_ids[i]] >= 0)
	groups[uid_to_group_2D[face_ids[i]]]->Append(i+1);
      

    outfile << "// Object2D GroupID, #Faces <immediately followed by> FaceID List:\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n";
    for(int i=0; i<numObj2D; i++)
      {
	outfile << i << " " << groups[i]->Size() << "\n";
	for(int j=0; j<groups[i]->Size(); j++)
	  {
	    outfile << (*groups[i])[j];
	    if(ports.Contains(face_ids[(*groups[i])[j]-1]))
	      outfile << " P";
	    outfile << "\n";
	  }
      }
    outfile << endl;

    
    for(int i=0; i<groups.Size(); i++)
      groups[i]->SetSize(0);

    for(int i=0; i<edge_ids.Size(); i++)
      if(uid_to_group_1D[edge_ids[i]] >= 0)
	groups[uid_to_group_1D[edge_ids[i]]]->Append(i+1);



    outfile << "// Object1D GroupID, #Edges <immediately followed by> EdgeID List:\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n";
    for(int i=0; i<numObj1D; i++)
      {
	outfile << i << " " << groups[i]->Size() << "\n";
	for(int j=0; j<groups[i]->Size(); j++)
	  outfile << (*groups[i])[j] << "\n";
      }
    outfile << endl;

    
    for(int i=0; i<groups.Size(); i++)
      groups[i]->SetSize(0);
    for(PointIndex i=PointIndex::BASE; i<mesh.GetNP()+PointIndex::BASE; i++)
      {
	if(i-PointIndex::BASE < point_ids.Size())
	  {
	    if(uid_to_group_0D[point_ids[i]] >= 0)
	      groups[uid_to_group_0D[point_ids[i]]]->Append(i+1-PointIndex::BASE);
	  }
	else
	  groups[uid_to_group_0D[0]]->Append(i+1-PointIndex::BASE);
      }


    outfile << "// Object0D GroupID, #Nodes <immediately followed by> NodeID List:\n" \
	    << "// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^\n";
    for(int i=0; i<numObj0D; i++)
      {
	outfile << i << " " << groups[i]->Size() << "\n";
	for(int j=0; j<groups[i]->Size(); j++)
	  outfile << (*groups[i])[j] << "\n";
      }
    outfile << endl;

    for(int i=0; i<groups.Size(); i++)
      delete groups[i];


    outfile.close();

    cout << ".tet export done" << endl;
  }
  void BaseMoveableMem :: Alloc (size_t s)
  {
    if (totalsize == 0)
      {
        size = s;
        //ptr = (char*) malloc(s);
        ptr = new char[s];

        if (!ptr)
          {
            cerr << "BaseynamicMem, cannot allocate " << s << " bytes" << endl;
            Print ();
            throw ("BaseDynamicMem::Alloc: out of memory");
          }

        return;
      }


    used += s - size;

    size_t r = s % 8;
    if (r) s += 8-r;
    if (prev)
      pos = prev->pos + prev->size;
    else
      pos = 0;
    size = s;
  
    if (next)
      {
        NgLock lock(mem_mutex);
        lock.Lock();
        try
          {
            next->MoveTo (pos+size);
          }
        catch (NgException e)
          {
            lock.UnLock();
            throw NgException ("MoveableMem overflow");
          }
        lock.UnLock();
      }

    if (size)
      {
        if (!largeblock)
          {
            cout << "moveable memory: allocate large block of "
                 << totalsize / 1048576  << " MB" << endl; 
            // largeblock = new char[totalsize];
            // largeblock = (char*)malloc (totalsize);
            largeblock = new char[totalsize];
          }
        ptr = largeblock+pos;

        if (pos + size > totalsize)
          throw NgException ("MoveableMem overflow");
      }
    else
      ptr = 0;
  }
Beispiel #10
0
void Primitive :: Transform (Transformation<3> & trans)
{
  stringstream ost;
  ost << "Primitve::Transform not implemented for " << typeid(*this).name() << endl;
  throw NgException (ost.str());
}
Beispiel #11
0
Primitive * Primitive :: Copy () const
{
  stringstream ost;
  ost << "Primitve::Copy not implemented for " << typeid(*this).name() << endl;
  throw NgException (ost.str());
}
Beispiel #12
0
 void NetgenGeometry :: Save (string filename) const
 {
   throw NgException("Cannot save geometry - no geometry available");
 }
Beispiel #13
0
  // extern double teterrpow; 
  MESHING3_RESULT MeshVolume (MeshingParameters & mp, Mesh& mesh3d)
  {
     int oldne;
     int meshed;

     Array<INDEX_2> connectednodes;

     if (&mesh3d.LocalHFunction() == NULL) mesh3d.CalcLocalH(mp.grading);

     mesh3d.Compress();

     //  mesh3d.PrintMemInfo (cout);

     if (mp.checkoverlappingboundary)
        if (mesh3d.CheckOverlappingBoundary())
           throw NgException ("Stop meshing since boundary mesh is overlapping");

     int nonconsist = 0;
     for (int k = 1; k <= mesh3d.GetNDomains(); k++)
     {
        PrintMessage (3, "Check subdomain ", k, " / ", mesh3d.GetNDomains());

        mesh3d.FindOpenElements(k);

        /*
        bool res = mesh3d.CheckOverlappingBoundary();
        if (res)
        {
        PrintError ("Surface is overlapping !!");
        nonconsist = 1;
        }
        */

        bool res = (mesh3d.CheckConsistentBoundary() != 0);
        if (res)
        {
           PrintError ("Surface mesh not consistent");
           nonconsist = 1;
        }
     }

     if (nonconsist)
     {
        PrintError ("Stop meshing since surface mesh not consistent");
        throw NgException ("Stop meshing since surface mesh not consistent");
     }

     double globmaxh = mp.maxh;

     for (int k = 1; k <= mesh3d.GetNDomains(); k++)
       {
	 if (multithread.terminate)
           break;
	 
	 PrintMessage (2, "");
	 PrintMessage (1, "Meshing subdomain ", k, " of ", mesh3d.GetNDomains());
	 (*testout) << "Meshing subdomain " << k << endl;
	 
	 mp.maxh = min2 (globmaxh, mesh3d.MaxHDomain(k));
	 
	 mesh3d.CalcSurfacesOfNode();
	 mesh3d.FindOpenElements(k);
	 
	 if (!mesh3d.GetNOpenElements())
           continue;
	 
	 

	 Box<3> domain_bbox( Box<3>::EMPTY_BOX ); 
	 
	 for (SurfaceElementIndex sei = 0; sei < mesh3d.GetNSE(); sei++)
	   {
	     const Element2d & el = mesh3d[sei];
	     if (el.IsDeleted() ) continue;
	     
	     if (mesh3d.GetFaceDescriptor(el.GetIndex()).DomainIn() == k ||
		 mesh3d.GetFaceDescriptor(el.GetIndex()).DomainOut() == k)
	       
	       for (int j = 0; j < el.GetNP(); j++)
		 domain_bbox.Add (mesh3d[el[j]]);
	   }
	 domain_bbox.Increase (0.01 * domain_bbox.Diam());
	 
	
        for (int qstep = 1; qstep <= 3; qstep++)
	  {
	    // cout << "openquads = " << mesh3d.HasOpenQuads() << endl;
	    if (mesh3d.HasOpenQuads())
	      {
		string rulefile = ngdir;
		
		const char ** rulep = NULL;
		switch (qstep)
		  {
		  case 1:
		    rulefile += "/rules/prisms2.rls";
		    rulep = prismrules2;
		    break;
		  case 2: // connect pyramid to triangle
		    rulefile += "/rules/pyramids2.rls";
		    rulep = pyramidrules2;
		    break;
		  case 3: // connect to vis-a-vis point
		    rulefile += "/rules/pyramids.rls";
		    rulep = pyramidrules;
		    break;
		  }
		
		//		Meshing3 meshing(rulefile);
		Meshing3 meshing(rulep); 
		
		MeshingParameters mpquad = mp;
		
		mpquad.giveuptol = 15;
		mpquad.baseelnp = 4;
		mpquad.starshapeclass = 1000;
		mpquad.check_impossible = qstep == 1;   // for prisms only (air domain in trafo)
		
		
		for (PointIndex pi = mesh3d.Points().Begin(); pi < mesh3d.Points().End(); pi++)
		  meshing.AddPoint (mesh3d[pi], pi);
		
		mesh3d.GetIdentifications().GetPairs (0, connectednodes);
		for (int i = 1; i <= connectednodes.Size(); i++)
		  meshing.AddConnectedPair (connectednodes.Get(i));
		
		for (int i = 1; i <= mesh3d.GetNOpenElements(); i++)
		  {
		    Element2d hel = mesh3d.OpenElement(i);
		    meshing.AddBoundaryElement (hel);
		  }
		
		oldne = mesh3d.GetNE();
		
		meshing.GenerateMesh (mesh3d, mpquad);
		
		for (int i = oldne + 1; i <= mesh3d.GetNE(); i++)
		  mesh3d.VolumeElement(i).SetIndex (k);
		
		(*testout) 
		  << "mesh has " << mesh3d.GetNE() << " prism/pyramid elements" << endl;
		
		mesh3d.FindOpenElements(k);
	      }
	  }
	

        if (mesh3d.HasOpenQuads())
        {
           PrintSysError ("mesh has still open quads");
           throw NgException ("Stop meshing since too many attempts");
           // return MESHING3_GIVEUP;
        }


        if (mp.delaunay && mesh3d.GetNOpenElements())
        {
           Meshing3 meshing((const char**)NULL);

           mesh3d.FindOpenElements(k);

           /*
           for (PointIndex pi = mesh3d.Points().Begin(); pi < mesh3d.Points().End(); pi++)
              meshing.AddPoint (mesh3d[pi], pi);
           */
           for (PointIndex pi : mesh3d.Points().Range())
              meshing.AddPoint (mesh3d[pi], pi);

           for (int i = 1; i <= mesh3d.GetNOpenElements(); i++)
              meshing.AddBoundaryElement (mesh3d.OpenElement(i));

           oldne = mesh3d.GetNE();

           meshing.Delaunay (mesh3d, k, mp);

           for (int i = oldne + 1; i <= mesh3d.GetNE(); i++)
              mesh3d.VolumeElement(i).SetIndex (k);

           PrintMessage (3, mesh3d.GetNP(), " points, ",
              mesh3d.GetNE(), " elements");
        }


        int cntsteps = 0;
        if (mesh3d.GetNOpenElements())
           do
           {
              if (multithread.terminate)
                 break;

              mesh3d.FindOpenElements(k);
              PrintMessage (5, mesh3d.GetNOpenElements(), " open faces");
              cntsteps++;

              if (cntsteps > mp.maxoutersteps) 
                 throw NgException ("Stop meshing since too many attempts");

              string rulefile = ngdir + "/tetra.rls";
              PrintMessage (1, "start tetmeshing");

              //	  Meshing3 meshing(rulefile);
              Meshing3 meshing(tetrules);

              Array<int, PointIndex::BASE> glob2loc(mesh3d.GetNP());
              glob2loc = -1;

              for (PointIndex pi = mesh3d.Points().Begin(); pi < mesh3d.Points().End(); pi++)
                if (domain_bbox.IsIn (mesh3d[pi]))
                    glob2loc[pi] = 
                    meshing.AddPoint (mesh3d[pi], pi);

              for (int i = 1; i <= mesh3d.GetNOpenElements(); i++)
              {
                 Element2d hel = mesh3d.OpenElement(i);
                 for (int j = 0; j < hel.GetNP(); j++)
                    hel[j] = glob2loc[hel[j]];
                 meshing.AddBoundaryElement (hel);
                 // meshing.AddBoundaryElement (mesh3d.OpenElement(i));
              }

              oldne = mesh3d.GetNE();

              mp.giveuptol = 15 + 10 * cntsteps; 
              mp.sloppy = 5;
              meshing.GenerateMesh (mesh3d, mp);

              for (ElementIndex ei = oldne; ei < mesh3d.GetNE(); ei++)
                 mesh3d[ei].SetIndex (k);


              mesh3d.CalcSurfacesOfNode();
              mesh3d.FindOpenElements(k);

              // teterrpow = 2;
              if (mesh3d.GetNOpenElements() != 0)
              {
                 meshed = 0;
                 PrintMessage (5, mesh3d.GetNOpenElements(), " open faces found");

                 MeshOptimize3d optmesh(mp);

                 const char * optstr = "mcmstmcmstmcmstmcm";
                 for (size_t j = 1; j <= strlen(optstr); j++)
                 {
                    mesh3d.CalcSurfacesOfNode();
                    mesh3d.FreeOpenElementsEnvironment(2);
                    mesh3d.CalcSurfacesOfNode();

                    switch (optstr[j-1])
                    {
                    case 'c': optmesh.CombineImprove(mesh3d, OPT_REST); break;
                    case 'd': optmesh.SplitImprove(mesh3d, OPT_REST); break;
                    case 's': optmesh.SwapImprove(mesh3d, OPT_REST); break;
                    case 't': optmesh.SwapImprove2(mesh3d, OPT_REST); break;
                    case 'm': mesh3d.ImproveMesh(mp, OPT_REST); break;
                    }	  

                 }

                 mesh3d.FindOpenElements(k);	      
                 PrintMessage (3, "Call remove problem");
                 RemoveProblem (mesh3d, k);
                 mesh3d.FindOpenElements(k);
              }
              else
              {
                 meshed = 1;
                 PrintMessage (1, "Success !");
              }
           }
           while (!meshed);

           PrintMessage (1, mesh3d.GetNP(), " points, ",
              mesh3d.GetNE(), " elements");
     }

     mp.maxh = globmaxh;

     MeshQuality3d (mesh3d);

     return MESHING3_OK;
  }  
Beispiel #14
0
void MeshOptimize2d :: ImproveVolumeMesh (Mesh & mesh)
{

    if (!faceindex)
    {
        PrintMessage (3, "Smoothing");

        for (faceindex = 1; faceindex <= mesh.GetNFD(); faceindex++)
        {
            ImproveVolumeMesh (mesh);
            if (multithread.terminate)
                throw NgException ("Meshing stopped");
        }
        faceindex = 0;
        return;
    }



    static int timer = NgProfiler::CreateTimer ("MeshSmoothing 2D");
    NgProfiler::RegionTimer reg (timer);



    CheckMeshApproximation (mesh);

    int i, j, k;
    SurfaceElementIndex sei;

    Array<SurfaceElementIndex> seia;
    mesh.GetSurfaceElementsOfFace (faceindex, seia);

    bool mixed = 0;
    for (i = 0; i < seia.Size(); i++)
        if (mesh[seia[i]].GetNP() != 3)
        {
            mixed = 1;
            break;
        }


    int loci;
    double fact;
    bool moveisok;

    PointGeomInfo ngi;
    Point<3> origp;

    Vector x(3);

    Array<MeshPoint, PointIndex::BASE> savepoints(mesh.GetNP());

    Array<int, PointIndex::BASE> nelementsonpoint(mesh.GetNP());
    nelementsonpoint = 0;

    for (i = 0; i < seia.Size(); i++)
    {
        const Element2d & el = mesh[seia[i]];
        for (j = 0; j < el.GetNP(); j++)
            nelementsonpoint[el[j]]++;
    }


    TABLE<SurfaceElementIndex,PointIndex::BASE> elementsonpoint(nelementsonpoint);
    for (i = 0; i < seia.Size(); i++)
    {
        const Element2d & el = mesh[seia[i]];
        for (j = 0; j < el.GetNP(); j++)
            elementsonpoint.Add (el[j], seia[i]);
    }


    JacobianPointFunction pf(mesh.Points(),mesh.VolumeElements());



//     Opti2SurfaceMinFunction surfminf(mesh);
//     Opti2EdgeMinFunction edgeminf(mesh);
//     Opti2SurfaceMinFunctionJacobian surfminfj(mesh);

    OptiParameters par;
    par.maxit_linsearch = 8;
    par.maxit_bfgs = 5;

    int np = mesh.GetNP();
    int ne = mesh.GetNE();

    BitArray badnodes(np);
    badnodes.Clear();

    for (i = 1; i <= ne; i++)
    {
        const Element & el = mesh.VolumeElement(i);
        double bad = el.CalcJacobianBadness (mesh.Points());
        if (bad > 1)
            for (j = 1; j <= el.GetNP(); j++)
                badnodes.Set (el.PNum(j));
    }


    bool printeddot = 0;
    char plotchar = '.';
    int modplot = 1;
    if (mesh.GetNP() > 1000)
    {
        plotchar = '+';
        modplot = 10;
    }
    if (mesh.GetNP() > 10000)
    {
        plotchar = 'o';
        modplot = 100;
    }
    int cnt = 0;


    Array<SurfaceElementIndex> locelements(0);
    Array<int> locrots(0);

    for (PointIndex pi = mesh.Points().Begin(); pi < mesh.Points().End(); pi++)
    {
        if (mesh[pi].Type() != SURFACEPOINT)
            continue;

        if (multithread.terminate)
            throw NgException ("Meshing stopped");

        int surfi(-1);

        if(elementsonpoint[pi].Size() == 0)
            continue;

        Element2d & hel = mesh[elementsonpoint[pi][0]];

        if(hel.GetIndex() != faceindex)
            continue;

        cnt++;
        if (cnt % modplot == 0 && writestatus)
        {
            printeddot = 1;
            PrintDot (plotchar);
        }


        int hpi = 0;
        for (j = 1; j <= hel.GetNP(); j++)
            if (hel.PNum(j) == pi)
            {
                hpi = j;
                break;
            }
        PointGeomInfo gi1 = hel.GeomInfoPi(hpi);

        locelements.SetSize(0);
        locrots.SetSize (0);

        for (j = 0; j < elementsonpoint[pi].Size(); j++)
        {
            sei = elementsonpoint[pi][j];
            const Element2d & bel = mesh[sei];
            surfi = mesh.GetFaceDescriptor(bel.GetIndex()).SurfNr();

            locelements.Append (sei);

            for (k = 1; k <= bel.GetNP(); k++)
                if (bel.PNum(k) == pi)
                {
                    locrots.Append (k);
                    break;
                }
        }


        double lh = mesh.GetH(mesh.Point(pi));
        par.typx = lh;

        pf.SetPointIndex(pi);

        x = 0;
        bool pok = (pf.Func (x) < 1e10);

        if (pok)
        {
            BFGS (x, pf, par);

            origp = mesh[pi];
            loci = 1;
            fact = 1;
            moveisok = false;


            //optimizer loop (if whole distance is not possible, move only a bit!!!!)
            while (loci <= 5 && !moveisok)
            {
                loci ++;
                mesh[pi](0) = origp(0) + x(0)*fact;
                mesh[pi](1) = origp(1) + x(1)*fact;
                mesh[pi](2) = origp(2) + x(2)*fact;
                fact = fact/2.;


                //cout << "origp " << origp << " newp " << mesh[pi];

                ngi = gi1;
                moveisok = (ProjectPointGI (surfi, mesh[pi], ngi) != 0);

                //cout << " projected " << mesh[pi] << endl;

                // point lies on same chart in stlsurface

                if (moveisok)
                {
                    for (j = 0; j < locelements.Size(); j++)
                        mesh[locelements[j]].GeomInfoPi(locrots[j]) = ngi;

                    //cout << "moved " << origp << " to " << mesh[pi] << endl;
                }
                else
                {
                    mesh[pi] = origp;
                }

            }
        }
        else
        {
            cout << "el not ok (point " << pi << ": " << mesh[pi] << ")" << endl;
        }
    }

    if (printeddot)
        PrintDot ('\n');

    CheckMeshApproximation (mesh);
    mesh.SetNextTimeStamp();
}