1 /* --------------------------------------------------------------------  */
2 /*                          CALCULIX                                     */
3 /*                   - GRAPHICAL INTERFACE -                             */
4 /*                                                                       */
5 /*     A 3-dimensional pre- and post-processor for finite elements       */
6 /*              Copyright (C) 1996 Klaus Wittig                          */
7 /*                                                                       */
8 /*     This program is free software; you can redistribute it and/or     */
9 /*     modify it under the terms of the GNU General Public License as    */
10 /*     published by the Free Software Foundation; version 2 of           */
11 /*     the License.                                                      */
12 /*                                                                       */
13 /*     This program is distributed in the hope that it will be useful,   */
14 /*     but WITHOUT ANY WARRANTY; without even the implied warranty of    */
15 /*     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the      */
16 /*     GNU General Public License for more details.                      */
17 /*                                                                       */
18 /*     You should have received a copy of the GNU General Public License */
19 /*     along with this program; if not, write to the Free Software       */
20 /*     Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.         */
21 /* --------------------------------------------------------------------  */
22 
23 /* do do:  */
24 /* merge point, check all nurbs if that point is used  */
25 
26 #include <cgx.h>
27 #define TEST 0
28 
29 
30 extern double     gtol;
31 
32 extern Scale     scale[1];
33 extern Summen    anz[1];
34 extern Edges     *edge;
35 
36 extern Nodes     *node;
37 extern Elements  *e_enqire;
38 extern Datasets *lcase;
39 extern Faces     *face;
40 
41 extern Alias     *alias;
42 extern Sets      *set;
43 extern Points    *point;
44 extern Lines     *line;
45 extern Lcmb      *lcmb;
46 extern Gsur      *surf;
47 extern Gbod      *body;
48 extern Nurbl     *nurbl;
49 extern Nurbs     *nurbs;
50 extern SumGeo    anzGeo[1];
51 
52 extern char  printFlag;                     /* printf 1:on 0:off */
53 extern char  addDispFlag;                    /* 0: original node-coordinates, 1: node-coordinates+displacements */
54 
55 
56 
57 /* additional entities from setFunktions */
58 
59 
60 
61 
62 /*---------------------------------------------------------------------*/
63 /* returns a field with sum and indexes of the entities inside a space */
64 /* of a radius of tol around mrgnod                                    */
65 /* type:  'n' node, 'p' point                                          */
66 /* rsort: radius and index of the entity, sorted by the radius         */
67 /* indx:  index of the entity for which the close ones are searched    */
68 /* anz_n: number of entities (indx is a member of rsort )              */
69 /* gtol:  geometric tolerance for which the search is done             */
70 /*---------------------------------------------------------------------*/
findCloseEntity(char type,Nodes * node,Points * point,Rsort * rsort,int indx,int anz_n,double local_gtol)71 int *findCloseEntity( char type, Nodes *node, Points *point, Rsort *rsort, int indx, int anz_n, double local_gtol)
72 {
73   int *nodes=NULL;
74   int n,m,i, n_low, n_high;
75   double rref, dx=0, dy=0, dz=0;
76 
77   /* find all nodes (nod_dr) which are in the intervall of rref-tol to rref+tol */
78   if( (nodes=(int *)malloc((1)*sizeof(int) ) )==NULL)
79   { printf(" ERROR: malloc failure\n"); return(NULL); }
80   nodes[0]=0;
81 
82   /* search the intersection on the lower bound */
83   rref=rsort[indx].r-local_gtol;
84   n_low=0;
85   n=anz_n;
86   do
87   {
88     m=(n+n_low)*.5;
89     if(rref >= rsort[m].r ) n_low=m;
90     else n=m;
91   }while((n-n_low) != 1);
92 
93   /* search the intersection on the upper bound */
94   rref=rsort[indx].r+local_gtol;
95   n_high=0;
96   n=anz_n;
97   do
98   {
99     m=(n+n_high)*.5;
100     if(rref >= rsort[m].r ) n_high=m;
101     else n=m;
102   }while((n-n_high) != 1);
103 
104 #if TEST
105   printf(" dn: %d n_low[%d]:%d n_high[%d]:%d\n", n_high-n_low, n_low,rsort[n_low].i, n_high,rsort[n_high].i);
106 #endif
107 
108   /* check all nodes inbetween */
109   for(i=n_low; i<=n_high; i++)
110   {
111     if((i!=indx)&&(rsort[i].i>-1))
112     {
113       /* is the entity relative to the master-entity in the tolerance-sphere? */
114       if(type=='n')
115       {
116         dx=node[rsort[indx].i].nx-node[rsort[i].i].nx; if(dx>local_gtol) goto nexti;
117         dy=node[rsort[indx].i].ny-node[rsort[i].i].ny; if(dy>local_gtol) goto nexti;
118         dz=node[rsort[indx].i].nz-node[rsort[i].i].nz; if(dz>local_gtol) goto nexti;
119       }
120       else if(type=='p')
121       {
122         dx=point[rsort[indx].i].px-point[rsort[i].i].px; if(dx>local_gtol) goto nexti;
123         dy=point[rsort[indx].i].py-point[rsort[i].i].py; if(dy>local_gtol) goto nexti;
124         dz=point[rsort[indx].i].pz-point[rsort[i].i].pz; if(dz>local_gtol) goto nexti;
125       }
126       if((dx*dx+dy*dy+dz*dz)<=local_gtol*local_gtol)
127       {
128         if( (nodes=(int *)realloc((int *)nodes, (nodes[0]+2)*sizeof(int) ) )==NULL)
129         { printf(" ERROR: realloc failure\n"); return(NULL); }
130         nodes[0]++;
131         nodes[nodes[0]]=i;
132       }
133 #if TEST
134       if(nodes[0])
135       {
136         printf(" entity:%d r:%lf \n", rsort[i].i, (dx*dx+dy*dy+dz*dz) );
137       }
138 #endif
139     }
140     nexti:;
141   }
142   return(nodes);
143 }
144 
145 /*------------------------------------------------------------------*/
146 /* merge node deletes that node if near one exists                */
147 /*------------------------------------------------------------------*/
148 
mergeNode(int * arrayIndx,int anz_n,int indx,int lock,int * nodes,int ** nod2elem)149 int mergeNode( int *arrayIndx, int anz_n, int indx, int lock, int *nodes, int **nod2elem )
150 {
151   int i, j, k,n,l,e;
152   int  ipuf, nod2, nod;
153 
154   nod=arrayIndx[indx];
155   if(( nod<1 )||( node[nod].pflag)) return(-1);
156 
157   for (i=1; i<=nodes[0]; i++)
158   {
159     if(nodes[i]<0) continue;
160     if (arrayIndx[nodes[i]] > 0)
161     {
162       nod2=arrayIndx[nodes[i]];
163       if(nod2<=0) continue;
164 
165       /* check if nod exists in related elements */
166       for (e=1; e<=nod2elem[nod2][0]; e++)
167       {
168         j=nod2elem[nod2][e];
169         // check if the elem still exists
170         if(e_enqire[j].type == 0) continue;
171         else if (e_enqire[j].type == 1) ipuf = 8;  /* CHEXA8  */
172         else if (e_enqire[j].type == 2) ipuf = 6; /* PENTA6 */
173         else if (e_enqire[j].type == 3) ipuf = 4;  /* CTET4   */
174         else if (e_enqire[j].type == 4) ipuf = 20; /* CHEXA20 */
175         else if (e_enqire[j].type == 5) ipuf = 15; /* PENTA15 */
176         else if (e_enqire[j].type == 6) ipuf = 10; /* CTET10  */
177         else if (e_enqire[j].type == 7) ipuf = 3;  /* CTRI3   */
178         else if (e_enqire[j].type == 8) ipuf = 6;  /* CTRI6   */
179         else if (e_enqire[j].type == 9) ipuf = 4;  /* CQUAD4  */
180         else if (e_enqire[j].type == 10) ipuf = 8;  /* CQUAD8  */
181         else if (e_enqire[j].type == 11) ipuf = 2; /* CBEAM2   */
182         else if (e_enqire[j].type == 12) ipuf = 3; /* CBEAM3   */
183         else ipuf=0;
184         for (k=0; k<ipuf; k++)
185         {
186           if( e_enqire[j].nod[k] == nod2 )
187           {
188             if (lock) for (l=0; l<ipuf; l++)
189             {
190               if( e_enqire[j].nod[l] == nod )
191               {
192                 printf(" node %d can not replace %d in element %d without collapsing it\n", nod, nod2, j );
193                 goto nextNode;
194               }
195             }
196             /* replace node */
197             if(printFlag) printf(" %d replace %d with %d in elem:%d type:%d\n", j, nod2, nod, j, e_enqire[j].type );
198             /* change node in elems if that node is used  */
199             e_enqire[j].nod[k] =nod;
200 
201             /* delete all elements which are defined by unique nodes */
202             if(!lock)
203             {
204               for(n=1; n<ipuf; n++) if(e_enqire[j].nod[n]!=e_enqire[e_enqire[i].nr].nod[0]) break;
205               if(n==ipuf)
206               { delElem(1,&j); j--; }
207 
208               /* if two nodes of a tr3 are unique */
209               if (e_enqire[j].type == 7) if((e_enqire[j].nod[0]==e_enqire[j].nod[1])||(e_enqire[j].nod[0]==e_enqire[j].nod[2])||(e_enqire[j].nod[1]==e_enqire[j].nod[2]))
210               { delElem(1,&j); j--; }
211 	    }
212           }
213         }
214       }
215 
216       /* check if nod exists in related entities */
217       for (j=0; j<anzGeo->p; j++) for (k=0; k<point[j].nn; k++) if( point[j].nod[k] == nod2 ) point[j].nod[k] =nod;
218       for (j=0; j<anzGeo->l; j++) for (k=0; k<line[j].nn; k++) if( line[j].nod[k] == nod2 ) line[j].nod[k] =nod;
219       for (j=0; j<anzGeo->s; j++) for (k=0; k<surf[j].nn; k++) if( surf[j].nod[k] == nod2 ) surf[j].nod[k] =nod;
220       for (j=0; j<anzGeo->b; j++) for (k=0; k<body[j].nn; k++) if( body[j].nod[k] == nod2 ) body[j].nod[k] =nod;
221       for(j=1; j<anz->sets; j++)
222       {
223         if(!set[j].type)
224         {
225           if((set[j].name!=(char *)NULL)&&( getIndex(&set[j].node,set[j].anz_n,nod2) >-1)) seta(j,"n",nod);
226 	}
227       }
228 
229       /* check if nod exists in faces */
230       /* -> not necessary because all faces will be destroyed and newly created after merge */
231 
232       delNod( 1, &nod2 );
233       arrayIndx[nodes[i]]=0;
234     }
235     nextNode:;
236   }
237 
238   return(1);
239 }
240 
241 /*------------------------------------------------------------------*/
242 /* merge elements                                                   */
243 /*------------------------------------------------------------------*/
244 
245 /* wird von qsort aufgerufen, vergleicht die ersten drei Zahlen von Integer-Feldern */
246 typedef struct {
247     int nr;
248     int nod[20];
249 }Esort;
250 
compareNodesOfElems(Esort * a,Esort * b)251 int compareNodesOfElems(Esort *a, Esort *b)
252 {
253   int i;
254 
255   /*
256     printf("a:");
257   for(i=19; i>=0; i--) printf(" %d", a[0].nod[i]);
258     printf("\n");
259     printf("b:");
260   for(i=19; i>=0; i--) printf(" %d", b[0].nod[i]);
261     printf("\n");
262   */
263 
264   for(i=19; i>=0; i--)
265   {
266     if(a[0].nod[i]>b[0].nod[i]) { return(1); }
267     if(a[0].nod[i]<b[0].nod[i]) { return(-1); }
268   }
269   return(0);
270 }
271 
mergeElem(int setNr)272 void mergeElem( int setNr )
273 {
274   int i, j,jj, k, indx, ipuf,n;
275   int sum_delem=0, *delem=NULL;
276 
277   int sum_e=0;
278   Esort *esort=NULL;
279 
280   /* sort all nodes of all elems in setNr */
281   for (i=0; i<set[setNr].anz_e; i++)
282   {
283     indx=set[setNr].elem[i];
284     if (e_enqire[indx].type == 1) ipuf = 8;  /* HEXA8 */
285     else if (e_enqire[indx].type == 2) ipuf = 6;  /* PENTA6 */
286     else if (e_enqire[indx].type == 3) ipuf = 4;  /* TET4 */
287     else if (e_enqire[indx].type == 4) ipuf = 20; /* HEX20 */
288     else if (e_enqire[indx].type == 5) ipuf = 15; /* PENTA15 */
289     else if (e_enqire[indx].type == 6) ipuf = 10; /* TET10 */
290     else if (e_enqire[indx].type == 7) ipuf = 3;  /* TRI3  */
291     else if (e_enqire[indx].type == 8) ipuf = 6;  /* TRI6  */
292     else if (e_enqire[indx].type == 9) ipuf = 4;  /* QUAD4 */
293     else if (e_enqire[indx].type == 10) ipuf = 8; /* QUAD8 */
294     else if (e_enqire[indx].type == 11) ipuf = 2; /* BEAM */
295     else if (e_enqire[indx].type == 12) ipuf = 3; /* BEAM3 */
296     else ipuf=0;
297 
298     if((esort=(Esort *)realloc((Esort *)esort, (sum_e+1)*sizeof(Esort)))==NULL)
299     { printf("\nERROR: realloc failure in mergElem()\n\n"); return; }
300     esort[i].nr=indx;
301     for (n=0; n<ipuf; n++) esort[sum_e].nod[n]=e_enqire[indx].nod[n];
302     for (; n<20; n++) esort[sum_e].nod[n]=0;
303 
304     /* 1st node has lowest number */
305     qsort( &esort[sum_e].nod[0], 20, sizeof(int), (void *)compareInt );
306     sum_e++;
307   }
308 
309   /* sort all sorted elems */
310   qsort( esort, sum_e, sizeof(Esort), (void *)compareNodesOfElems );
311 
312   /* delete all elems which are equal to the last elem */
313   j=0; i=1; while(i<set[setNr].anz_e)
314   {
315     n=compareNodesOfElems(&esort[j],&esort[i]);
316     if(n==0)
317     {
318       /* Check if esort[i].nr exists in related entities and replace it by esort[j] */
319       for (jj=0; jj<anzGeo->l; jj++) for (k=0; k<line[jj].ne; k++) if( line[jj].elem[k] == esort[i].nr ) line[jj].elem[k] =esort[j].nr;
320       for (jj=0; jj<anzGeo->s; jj++) for (k=0; k<surf[jj].ne; k++) if( surf[jj].elem[k] == esort[i].nr ) surf[jj].elem[k] =esort[j].nr;
321       for (jj=0; jj<anzGeo->b; jj++) for (k=0; k<body[jj].ne; k++) if( body[jj].elem[k] == esort[i].nr ) body[jj].elem[k] =esort[j].nr;
322       for(jj=1; jj<anz->sets; jj++)
323       {
324         if(!set[jj].type)
325         {
326           if((set[jj].name!=(char *)NULL)&&( getIndex(&set[jj].elem,set[jj].anz_e, esort[i].nr) >-1)) seta(jj,"e",esort[j].nr);
327 	}
328       }
329 
330       if((delem=(int *)realloc((int *)delem, (sum_delem+1)*sizeof(int)))==NULL)
331       { printf("\nERROR: realloc failure in mergElem()\n\n"); return; }
332       delem[sum_delem]=esort[i].nr;
333       sum_delem++;
334     }
335     else j=i;
336     i++;
337   };
338 
339 
340   /* delete element */
341   delElem( sum_delem, delem);
342   printf("%d elems deleted\n", sum_delem);
343 
344   if(delem) free(delem);
345   if(esort) free(esort);
346   return;
347 }
348 
349 /*------------------------------------------------------------------*/
350 /* merge point deletes all matching points                          */
351 /*------------------------------------------------------------------*/
352 
mergePnt(Rsort * rsort,int anz_p,int indx,double local_gtol,int lock)353 int mergePnt( Rsort *rsort, int anz_p, int indx, double local_gtol, int lock )
354 {
355   int i, j, k, n;
356   int *pnts, pnt, pnt2;
357 
358   pnt=rsort[indx].i;
359   if ( pnt<0 ) return(-1);
360   if ( anz_p<2 ) return(-1);
361 
362 #if TEST
363   printf("  pnt[%d]:%s indx:%d anz_p:%d\n", pnt, point[pnt].name, indx, anz_p);
364 #endif
365 
366   pnts= findCloseEntity('p', 0, point, rsort, indx, anz_p, local_gtol);
367 
368   for (i=1; i<=pnts[0]; i++)
369   {
370     if (rsort[pnts[i]].i > -1)
371     {
372       pnt2=rsort[pnts[i]].i;
373 
374       /* check all lines and lcmbs if that point is used  */
375       if(lock) for (j=0; j<anzGeo->l; j++)
376       {
377         if( line[j].name != (char *)NULL )
378         {
379           if (line[j].p1 == pnt2)
380           {
381             if (line[j].p2==pnt)
382             {
383               printf(" P1:%s not changed to %s, would be equal to P2:%s in line:%s\n",
384                 point[line[j].p1].name, point[pnt].name, point[line[j].p2].name, line[j].name);
385               return(-1);
386             }
387           }
388           if (line[j].p2 == pnt2)
389           {
390             if (line[j].p1==pnt)
391             {
392               printf(" P2:%s not changed to %s, would be equal to P1:%s in line:%s\n",
393                 point[line[j].p2].name, point[pnt].name, point[line[j].p1].name, line[j].name);
394               return(-1);
395             }
396           }
397           if (line[j].typ=='a')
398 	  {
399             if ((line[j].p1 == pnt2)||(line[j].p2 == pnt2))
400             {
401               if (line[j].trk==pnt)
402               {
403                 printf(" P:%s not changed to %s, would be equal to cp:%s in line:%s\n",
404                   point[pnt2].name, point[pnt].name, point[line[j].trk].name, line[j].name);
405                 return(-1);
406               }
407             }
408             if ((line[j].p1 == pnt)||(line[j].p2 == pnt))
409             {
410               if (line[j].trk==pnt2)
411               {
412                 printf(" cp:%s not changed to %s, would be equal to P:%s in line:%s\n",
413                    point[pnt2].name, point[pnt].name, point[pnt].name, line[j].name);
414                 return(-1);
415               }
416             }
417           }
418         }
419       }
420 
421       if(lock) for (j=0; j<anzGeo->c; j++)
422       {
423         if( lcmb[j].name != (char *)NULL )
424         {
425           if (lcmb[j].p1 == pnt2)
426           {
427             if (lcmb[j].p2==pnt)
428             {
429               printf(" P1:%s not changed to %s, would be equal to P2:%s in lcmb:%s\n",
430                 point[lcmb[j].p1].name, point[pnt].name, point[lcmb[j].p2].name, lcmb[j].name);
431               return(-1);
432             }
433           }
434           if (lcmb[j].p2 == pnt2)
435           {
436             if (lcmb[j].p1==pnt)
437             {
438               printf(" P2:%s not changed to %s, would be equal to P1:%s in lcmb:%s\n",
439                 point[lcmb[j].p2].name, point[pnt].name, point[lcmb[j].p1].name, lcmb[j].name);
440               return(-1);
441             }
442           }
443         }
444       }
445 
446 
447       /* loesche point */
448       printf(" replace %s with %s\n", point[pnt2].name, point[pnt].name );
449 
450       /* change point in  lines if that point is used  */
451       for (j=0; j<anzGeo->l; j++)
452       {
453         if( line[j].name != (char *)NULL)
454         {
455           if (line[j].p1 == pnt2) line[j].p1=pnt;
456           if (line[j].p2 == pnt2) line[j].p2=pnt;
457           if ((line[j].trk == pnt2)&&(line[j].typ == 'a')) line[j].trk=pnt;
458         }
459       }
460       /* change point in  lcmbs if that point is used  */
461       for (j=0; j<anzGeo->c; j++)
462       {
463         if( lcmb[j].name != (char *)NULL)
464         {
465           if (lcmb[j].p1 == pnt2) lcmb[j].p1=pnt;
466           if (lcmb[j].p2 == pnt2) lcmb[j].p2=pnt;
467         }
468       }
469       /* change point in nurbs lines if that point is used  */
470       for (j=0; j<anzGeo->nurl; j++)
471       {
472         if( nurbl[j].name != (char *)NULL)
473         {
474           for (n=0; n<nurbl[j].u_npnt; n++)
475             if ( nurbl[j].ctlpnt[n]== pnt2) nurbl[j].ctlpnt[n]=pnt;
476         }
477       }
478       /* change point in nurbs-surfs if that point is used  */
479       for (j=0; j<anzGeo->nurs; j++)
480       {
481         if( nurbs[j].name != (char *)NULL)
482         {
483           for (n=0; n<nurbs[j].u_npnt; n++)
484             for (k=0; k<nurbs[j].v_npnt; k++)
485               if ( nurbs[j].ctlpnt[n][k]== pnt2) nurbs[j].ctlpnt[n][k]=pnt;
486         }
487       }
488       /* add point in set if pnt2 is a member and if it is a spline-set */
489       for (j=0; j<anz->sets; j++)
490       {
491         if( set[j].name != (char *)NULL)
492 	{
493           if( set[j].type==1)
494           {
495             for (k=0; k<set[j].anz_p; k++)
496             {
497               if (set[j].pnt[k] == pnt2)
498               {
499                 /* point to delete is a member of this set */
500                 /* check if the replace point is already a member */
501                 /* if not replace the point index by the replace point */
502                 for (n=0; n<set[j].anz_p; n++)
503                 {
504                   if( set[j].pnt[n] == pnt ) goto nextset;
505                 }
506                 set[j].pnt[k]=pnt; goto nextset;
507               }
508             }
509           }
510           else
511           {
512             if( getIndex(&set[j].pnt,set[j].anz_p,pnt2) >-1) seta(j,"p",pnt);
513           }
514           nextset:;
515         }
516       }
517       delPnt( 1, &pnt2 ); rsort[pnts[i]].i=-1;
518     }
519   }
520   return(1);
521 }
522 
523 /*------------------------------------------------------------------*/
524 /* merge line deletes all matching lines                           */
525 /*------------------------------------------------------------------*/
526 
mergeLine(int indx,int lock)527 int mergeLine( int indx, int lock)
528 {
529   int i, j, k, l, delflag, equalFlag;
530 
531   if( line[indx].name == (char *)NULL ) return(-1);
532 
533   /* check if a line has the same end-points */
534 
535   delflag=1;
536   for (i=0; i<anzGeo->l; i++)
537   {
538    if(( line[i].name != (char *)NULL ) && (i!=indx))
539    {
540     /* is the line principaly deletable? */
541 
542     if (line[indx].typ == line[i].typ)
543     if ( ((line[indx].p1 == line[i].p1)&&(line[indx].p2 == line[i].p2)) ||
544          ((line[indx].p1 == line[i].p2)&&(line[indx].p2 == line[i].p1)) )
545     {
546       equalFlag=0;
547       if(line[indx].trk != line[i].trk)
548       {
549         /* check if the trks are identicaly formulatet */
550         if ((line[indx].typ=='s')&&(line[i].typ=='s'))
551         {
552           if(set[line[indx].trk].anz_p==set[line[i].trk].anz_p)
553           {
554             equalFlag=1;
555             /* check forward */
556             for(j=0; j<set[line[indx].trk].anz_p; j++)
557             {   if(set[line[indx].trk].pnt[j]!=set[line[i].trk].pnt[j]) equalFlag=0; }
558             /* check backward */
559             if(!equalFlag)
560             {
561               equalFlag=1;
562               for(j=0; j<set[line[indx].trk].anz_p; j++)
563               {   if(set[line[indx].trk].pnt[j]!=set[line[i].trk].pnt[set[line[indx].trk].anz_p-1-j]) equalFlag=0; }
564             }
565           }
566         }
567       }
568       if(((line[indx].trk == line[i].trk)||(equalFlag))||(!lock))
569       {
570        /* check lcmbs for the line[indx]  */
571        for (j=0; j<anzGeo->c; j++)
572        {
573         if( lcmb[j].name != (char *)NULL) for (k=0; k<lcmb[j].nl; k++)
574         {
575           if (lcmb[j].l[k] == indx)
576           {
577             /* check if the line i would appear twice  */
578             for (l=0; l<lcmb[j].nl; l++)
579               if (lcmb[j].l[l]==i)
580               {
581                 errMsg (" can not merge %s, the line would appear twice in lcmb %d\n",
582 			line[indx].name, lcmb[j].l[l]);
583                 delflag=0;
584               }
585           }
586         }
587        }
588       }
589       else { printf(" can not merge line:%s with line:%s. they use different tracks\n",line[indx].name,line[i].name); delflag=0; }
590 
591       if (delflag)
592       {
593         /* check if a surf uses both identical lines line[indx] and line[i] and if yes remove i */
594         for (j=0; j<anzGeo->s; j++)
595         {
596           if( surf[j].name != (char *)NULL) for (k=0; k<surf[j].nl; k++)
597           {
598             if ((surf[j].l[k] == indx)&&(surf[j].typ[k]=='l')&&( surf[j].name != (char *)NULL))
599             {
600               /* check if the line i would appear twice  */
601               for (l=0; l<surf[j].nl; l++)
602               {
603                 if ((surf[j].typ[l]=='l')&&(surf[j].l[l]==i))
604                 {
605                   k=0;
606                   /* identical lines indx and i appear both in surface, delete i from the definition */
607                   for (l=0; l<surf[j].nl; l++)
608                   {
609                     if (surf[j].l[l] == i) continue;
610                     surf[j].l[k] = surf[j].l[l];
611                     surf[j].typ[k] = surf[j].typ[l];
612                     surf[j].o[k] = surf[j].o[l];
613                     k++;
614                   }
615                   surf[j].nl=k;
616                   goto let_surfs_collapse;
617                 }
618               }
619             }
620           }
621         }
622     let_surfs_collapse:;
623 
624       /* delete the "i" line */
625         printf(" %s replace with %s\n", line[i].name, line[indx].name );
626         /* check for equal oriented lines  */
627         if ((line[indx].p1 == line[i].p1)&&(line[indx].p2 == line[i].p2))
628         {
629           /* check lcmbs for the line[i]  */
630           for (j=0; j<anzGeo->c; j++)
631           {
632             if( lcmb[j].name != (char *)NULL) for (k=0; k<lcmb[j].nl; k++)
633             {
634               if (lcmb[j].l[k] == i) lcmb[j].l[k]=indx;
635             }
636           }
637           /* check surfs for the line[i]  */
638           for (j=0; j<anzGeo->s; j++)
639           {
640             if( surf[j].name != (char *)NULL) for (k=0; k<surf[j].nl; k++)
641             {
642               if ((surf[j].l[k] == i)&&(surf[j].typ[k]=='l')) surf[j].l[k]=indx;
643             }
644           }
645         }
646         /* line is reverse oriented */
647         else
648         {
649           /* check lcmbs for the line[i]  */
650           for (j=0; j<anzGeo->c; j++)
651           {
652             if( lcmb[j].name != (char *)NULL) for (k=0; k<lcmb[j].nl; k++)
653             {
654               if (lcmb[j].l[k] == i)
655               {
656                 lcmb[j].l[k]=indx;
657                 if(lcmb[j].o[k]=='+') lcmb[j].o[k]='-';
658                 else                 lcmb[j].o[k]='+';
659               }
660             }
661           }
662           /* check surfs for the line[i]  */
663           for (j=0; j<anzGeo->s; j++)
664           {
665             if( surf[j].name != (char *)NULL) for (k=0; k<surf[j].nl; k++)
666             {
667               if ((surf[j].l[k] == i)&&(surf[j].typ[k]=='l'))
668               {
669                 surf[j].l[k]=indx;
670                 if(surf[j].o[k]=='+') surf[j].o[k]='-';
671                 else                 surf[j].o[k]='+';
672               }
673             }
674           }
675         }
676         /* Check if line exists in related sets */
677         for (j=0; j<anz->sets; j++)
678         {
679           if(( set[j].name != (char *)NULL)&&( set[j].type==0))
680 	  {
681             if( getIndex(&set[j].line,set[j].anz_l,i) >-1) seta(j,"l",indx);
682           }
683         }
684 
685         delLine( 1, &i );
686       }
687     }
688    }
689   }
690   return(0);
691 }
692 
693 
694 /*------------------------------------------------------------------*/
695 /* merge lcmb deletes all matching lcmbs                            */
696 /*------------------------------------------------------------------*/
697 
mergeLcmb(int indx)698 int mergeLcmb( int indx )
699 {
700   int i, j, k, l;
701 
702   if( lcmb[indx].name == (char *)NULL) return(-1);
703 
704   /* check if a lcmb has the same composition */
705 
706   for (i=0; i<anzGeo->c; i++)
707   {
708     if(( lcmb[i].name != (char *)NULL) && (i!=indx) )
709     {
710       /* is the lcmb principaly deletable?  */
711 
712       if( ((lcmb[indx].p1 == lcmb[i].p1)&&(lcmb[indx].p2 == lcmb[i].p2)) ||
713          ((lcmb[indx].p1 == lcmb[i].p2)&&(lcmb[indx].p2 == lcmb[i].p1)) )
714       {
715         /* check surfs for the lcmb[i] and lcmb[indx] */
716         for (j=0; j<anzGeo->s; j++) if( surf[j].name != (char *)NULL)
717         {
718           for (k=0; k<surf[j].nl; k++)
719           {
720             if ((surf[j].l[k] == indx)&&(surf[j].typ[k]=='c'))
721             {
722               /* check if the lcmb indx would appear twice in a surf */
723               for (l=0; l<surf[j].nl; l++)
724                 if ((surf[j].typ[l]=='c')&&(surf[j].l[l]==i))
725                 {
726                   /* a surface would colapse, delete the surface */
727                   delSurf( 1, &j );
728                 }
729             }
730           }
731         }
732 
733         /* loesche lcmb */
734         printf(" %s replace with %s\n", lcmb[i].name, lcmb[indx].name );
735 
736         /* check surfs for the lcmb[i]  */
737         for (j=0; j<anzGeo->s; j++) if( surf[j].name != (char *)NULL )
738         {
739           for (k=0; k<surf[j].nl; k++)
740           {
741             if ((surf[j].l[k] == i)&&(surf[j].typ[k]=='c'))
742             {
743               if((lcmb[indx].p1 == lcmb[i].p1)&&(lcmb[indx].p2 == lcmb[i].p2))
744               {
745                 surf[j].l[k]=indx;
746               }
747               /* lcmb is reverse oriented */
748               else
749               {
750                 surf[j].l[k]=indx;
751                 if(surf[j].o[k]=='+') surf[j].o[k]='-';
752                 else                 surf[j].o[k]='+';
753               }
754             }
755           }
756         }
757         /* Check if lcmb exists in related sets */
758         for (j=0; j<anz->sets; j++)
759         {
760           if(( set[j].name != (char *)NULL)&&( set[j].type==0))
761 	  {
762             if( getIndex(&set[j].lcmb,set[j].anz_c,i) >-1) seta(j,"c",indx);
763           }
764         }
765 
766         delLcmb( 1, &i );
767         return(i);
768       }
769     }
770   }
771   return(-1);
772 }
773 
774 
775 /*------------------------------------------------------------------*/
776 /* merge surface                                                    */
777 /*------------------------------------------------------------------*/
778 
mergeSurf(int indx,double local_gtol)779 int mergeSurf( int indx, double local_gtol )
780 {
781   int i, j, k, n, l, rev;
782   int *found=NULL;
783   double r_indx, r_i;
784 
785   if( surf[indx].name == (char *)NULL) return(-1);
786 
787   if ( (found = (int *)realloc( (int *)found,(surf[indx].nl+1) * sizeof(int))) == NULL )
788   { printf("\n\n ERROR: realloc failed in mergeSurf()\n"); return(-1); }
789 
790   /* check if a surf has the same composition */
791 
792   for (i=0; i<anzGeo->s; i++)
793   {
794     /* is the surf principaly deletable?  */
795 
796     if( surf[i].name == (char *)NULL) goto nextsi;
797     if( i==indx ) goto nextsi;
798     if (surf[indx].nl!=surf[i].nl) goto nextsi;
799 
800     /* calculate the distance between the surfaces */
801     r_indx=sqrt(surf[indx].cx*surf[indx].cx+surf[indx].cy*surf[indx].cy+surf[indx].cz*surf[indx].cz);
802     r_i   =sqrt(surf[i].cx*surf[i].cx+surf[i].cy*surf[i].cy+surf[i].cz*surf[i].cz);
803     if ( abs(r_indx-r_i)> local_gtol) goto nextsi;
804 
805     /* compare the embedded indexes and types */
806     for (j=0; j<surf[indx].nl; j++)
807     {
808       found[j]=-1;
809       for (k=0; k<surf[i].nl; k++)
810       {
811         if( (surf[indx].l[j] == surf[i].l[k])&&(surf[indx].typ[j] == surf[i].typ[k]) ) found[j]=k;
812       }
813       if (found[j]==-1) goto nextsi;
814     }
815 
816     /* the surf can be merged */
817     /* check bodys for the surf[indx] and surf[i]  */
818     for (j=0; j<anzGeo->b; j++)
819     {
820       for (k=0; k<body[j].ns; k++)
821       {
822         if (( body[j].name != (char *)NULL )&&(body[j].s[k] == indx))
823         {
824            /* check if the surf indx would appear twice in a body */
825            for (l=0; l<body[j].ns; l++)
826              if (body[j].s[l]==i)
827              {
828                 /* a body would colapse, delete the body */
829                 delBody( 1, &j );
830              }
831         }
832       }
833     }
834 
835     /* loesche surf */
836     printf(" %s replace with %s\n", surf[i].name, surf[indx].name );
837 
838     /* check bodys for the surf[i]  */
839     for (j=0; j<anzGeo->b; j++) if( body[j].name != (char *)NULL )
840     {
841       for (k=0; k<body[j].ns; k++)
842       {
843         if (body[j].s[k] == i)
844         {
845           /* check the orientation of the entities in surf to replace */
846           rev=0;
847           for (n=0; n<surf[i].nl; n++)
848           {
849             if(surf[indx].o[n] != surf[i].o[found[n]]) rev=1;
850           }
851 	  if (!rev) /* the entities of surf i are oriented as in indx */
852           {
853             /* check the orientation of the surf to replace */
854             if (surf[indx].ori!=surf[i].ori)
855 	    {
856 	      if (body[j].o[k]=='+') body[j].o[k]='-';
857 	      else body[j].o[k]='+';
858 	    }
859             body[j].s[k]=indx;
860           }
861 	  if (rev) /* the entities of surf i are reverse oriented as in indx */
862           {
863             /* check the orientation of the surf to replace */
864             if (surf[indx].ori==surf[i].ori)
865 	    {
866 	      if (body[j].o[k]=='+') body[j].o[k]='-';
867 	      else body[j].o[k]='+';
868 	    }
869             body[j].s[k]=indx;
870           }
871         }
872       }
873     }
874     /* Check sets for the surf[i] */
875     for (j=0; j<anz->sets; j++)
876     {
877       if(( set[j].name != (char *)NULL)&&( set[j].type==0))
878       {
879         if( getIndex(&set[j].surf,set[j].anz_s,i) >-1) seta(j,"s",indx);
880       }
881     }
882 
883     delSurf( 1, &i );
884     free(found);
885     return(i);
886 
887  nextsi:;
888   }
889   free(found);
890   return(-1);
891 }
892 
893 
894 /*------------------------------------------------------------------*/
895 /* merge body                                                      */
896 /*------------------------------------------------------------------*/
897 
898 /*------------------------------------------------------------------*/
899 /* merge Nurbs surface                                             */
900 /*------------------------------------------------------------------*/
901 
902 
903 
pre_merge(char * record)904 void pre_merge( char *record )
905 {
906   int i,j,k,l,n,s;
907   int *indx=NULL;
908   int  sum, lock=1;
909   int  length, setNr;
910   char type[MAX_LINE_LENGTH], buffer[MAX_LINE_LENGTH], setname[MAX_LINE_LENGTH], alock[MAX_LINE_LENGTH];
911   double local_gtol=0.;
912   char  addDispFlagLocal=0;
913   Rsort *rsort=NULL;
914 
915   int     n_closest_nodes=0, ipuf=0, counter;
916   double *orig_x=NULL, *orig_y=NULL, *orig_z=NULL, *sort_x=NULL, *sort_y=NULL, *sort_z=NULL;
917   int *sort_nx=NULL, *sort_ny=NULL, *sort_nz=NULL, *near_node=NULL, *arrayIndx=NULL;
918   double dx=0., dy=0., dz=0.;
919   int **nod2elem=NULL;
920 
921   alock[0]=' ';
922   length= sscanf( record, "%s%s%s%s", type, setname, buffer, alock);
923 
924   if (length<2) return;
925 
926   operateAlias( setname, "se" );
927   setNr=getSetNr(setname);
928   if (setNr<0)
929   {
930     printf (" merge: set:%s does not exist\n", setname);
931     return;
932   }
933   printf(" merge ongoing, please wait\n");
934 
935   local_gtol=gtol;
936   if (length==3)
937   {
938     if(type[0]=='l')
939     { if ((buffer[0]=='n')||(buffer[0]=='N')) lock=0; }
940     else
941     {
942       if(!checkIfNumber(buffer)) { if ((buffer[0]=='n')||(buffer[0]=='N')) lock=0; }
943       else local_gtol=atof(buffer);
944     }
945   }
946   if (length==4)
947   {
948     local_gtol=atof(buffer);
949     if ((alock[0]=='n')||(alock[0]=='N')) lock=0;
950   }
951   /* since the coordinates are scaled: */
952   local_gtol/=scale->w;
953 
954   if (( type[0]=='n' )||( type[0]=='N' ))
955   {
956     if(addDispFlag==1)
957     {
958       addDispToCoordinates(node);
959       // remember to switch back
960       addDispFlagLocal=2;
961     }
962 
963     /* generate a table which relates elements to nodes */
964     if((nod2elem = (int **)malloc((anz->nmax+1) * sizeof(int *))) == NULL )
965       printf("ERROR: malloc failed\n\n");
966     for(i=0; i<=anz->nmax; i++)
967     {
968       if((nod2elem[i] = (int *)malloc((1) * sizeof(int))) == NULL )
969         printf("ERROR: malloc failed\n\n");
970       nod2elem[i][0]=0;
971     }
972 
973     for(i=0; i<anz->e; i++)
974     {
975       if (e_enqire[e_enqire[i].nr].type == 1) ipuf = 8;       /* HEXA8 */
976       else if (e_enqire[e_enqire[i].nr].type == 2) ipuf = 6;  /* PENTA6 */
977       else if (e_enqire[e_enqire[i].nr].type == 3) ipuf = 4;  /* TET4 */
978       else if (e_enqire[e_enqire[i].nr].type == 4) ipuf = 20; /* HEXA20 */
979       else if (e_enqire[e_enqire[i].nr].type == 5) ipuf = 15; /* PENTA15 */
980       else if (e_enqire[e_enqire[i].nr].type == 6) ipuf = 10; /* TET10 */
981       else if (e_enqire[e_enqire[i].nr].type == 7) ipuf = 3;  /* TRI3  */
982       else if (e_enqire[e_enqire[i].nr].type == 8) ipuf = 6;  /* TRI6  */
983       else if (e_enqire[e_enqire[i].nr].type == 9) ipuf = 4;  /* QUAD4 */
984       else if (e_enqire[e_enqire[i].nr].type == 10) ipuf = 8; /* QUAD8 */
985       else if (e_enqire[e_enqire[i].nr].type == 11) ipuf = 2; /* BEAM */
986       else if (e_enqire[e_enqire[i].nr].type == 12) ipuf = 3; /* BEAM3 */
987       else ipuf=0;
988       for (j=0; j<ipuf; j++)
989       {
990         n=e_enqire[e_enqire[i].nr].nod[j];
991 	if((nod2elem[n] = (int *)realloc((int *)nod2elem[n], (nod2elem[n][0]+2) * sizeof(int))) == NULL )
992           printf("ERROR: realloc failed\n\n");
993         nod2elem[n][0]++;
994         nod2elem[n][nod2elem[n][0]]=e_enqire[i].nr;
995       }
996     }
997 
998     /* sort the nodes */
999     /* stelle daten fuer near3d bereit */
1000     if ( (arrayIndx = (int *)malloc( (set[setNr].anz_n) * sizeof(int))) == NULL )
1001       printf("ERROR: malloc failed: Rsort\n\n" );
1002     if ( (near_node = (int *)malloc( (set[setNr].anz_n) * sizeof(int))) == NULL )
1003       printf("ERROR: malloc failed: Rsort\n\n" );
1004     if ( (rsort = (Rsort *)malloc( (set[setNr].anz_n) * sizeof(Rsort))) == NULL )
1005       printf("ERROR: malloc failed: Rsort\n\n" );
1006     if ( (orig_x = (double *)malloc( (set[setNr].anz_n) * sizeof(double))) == NULL )
1007       printf("ERROR: malloc failed in areampc\n\n" );
1008     if ( (orig_y = (double *)malloc( (set[setNr].anz_n) * sizeof(double))) == NULL )
1009       printf("ERROR: malloc failed in areampc\n\n" );
1010     if ( (orig_z = (double *)malloc( (set[setNr].anz_n) * sizeof(double))) == NULL )
1011       printf("ERROR: malloc failed in areampc\n\n" );
1012     if ( (sort_x = (double *)malloc( (set[setNr].anz_n) * sizeof(double))) == NULL )
1013       printf("ERROR: malloc failed in areampc\n\n" );
1014     if ( (sort_y = (double *)malloc( (set[setNr].anz_n) * sizeof(double))) == NULL )
1015       printf("ERROR: malloc failed in areampc\n\n" );
1016     if ( (sort_z = (double *)malloc( (set[setNr].anz_n) * sizeof(double))) == NULL )
1017       printf("ERROR: malloc failed in areampc\n\n" );
1018     if ( (sort_nx = (int *)malloc( (set[setNr].anz_n) * sizeof(int))) == NULL )
1019       printf("ERROR: malloc failed in areampc\n\n" );
1020     if ( (sort_ny = (int *)malloc( (set[setNr].anz_n) * sizeof(int))) == NULL )
1021       printf("ERROR: malloc failed in areampc\n\n" );
1022     if ( (sort_nz = (int *)malloc( (set[setNr].anz_n) * sizeof(int))) == NULL )
1023       printf("ERROR: malloc failed in areampc\n\n" );
1024 
1025     for(i=0; i<set[setNr].anz_n; i++)
1026     {
1027       arrayIndx[i]=set[setNr].node[i];
1028       rsort[i].r=orig_x[i]=node[set[setNr].node[i]].nx;
1029       rsort[i].i=i;
1030     }
1031     qsort( rsort, set[setNr].anz_n, sizeof(Rsort), (void *)compareRsort );
1032     for(i=0; i<set[setNr].anz_n; i++)
1033     {
1034       sort_x[i]=rsort[i].r;
1035       sort_nx[i]=rsort[i].i;
1036     }
1037     for(i=0; i<set[setNr].anz_n; i++)
1038     {
1039       rsort[i].r=orig_y[i]=node[set[setNr].node[i]].ny;
1040       rsort[i].i=i;
1041     }
1042     qsort( rsort, set[setNr].anz_n, sizeof(Rsort), (void *)compareRsort );
1043     for(i=0; i<set[setNr].anz_n; i++)
1044     {
1045       sort_y[i]=rsort[i].r;
1046       sort_ny[i]=rsort[i].i;
1047     }
1048     for(i=0; i<set[setNr].anz_n; i++)
1049     {
1050       rsort[i].r=orig_z[i]=node[set[setNr].node[i]].nz;
1051       rsort[i].i=i;
1052     }
1053     qsort( rsort, set[setNr].anz_n, sizeof(Rsort), (void *)compareRsort );
1054     for(i=0; i<set[setNr].anz_n; i++)
1055     {
1056       sort_z[i]=rsort[i].r;
1057       sort_nz[i]=rsort[i].i;
1058     }
1059 
1060     k=set[setNr].anz_n;
1061     counter=1e3;
1062     for(i=0; i<k; i++)
1063     {
1064       if(arrayIndx[i]==0) continue;
1065       n_closest_nodes=10;
1066       if(i>counter) { printf(" evaluated:%d from:%d\n", i, k); counter*=10; }
1067     repeatNodSearch:;
1068       if(n_closest_nodes>=k) n_closest_nodes=k-1;
1069       near3d(orig_x,orig_y,orig_z,sort_x,sort_y,sort_z,sort_nx,sort_ny,sort_nz, node[arrayIndx[i]].nx,node[arrayIndx[i]].ny,
1070         node[arrayIndx[i]].nz, k, &near_node[1], n_closest_nodes);
1071       /* store only nodes local_gtol away from node:arrayIndx[i] */
1072       /* do we have all nodes in range? */
1073       if(n_closest_nodes<k-1)
1074       {
1075         // arrayIndx[near_node[n_closest_nodes]]==nodeNr, if deleted set to 0 in mergNode()
1076         ipuf=n_closest_nodes;
1077         while(arrayIndx[near_node[ipuf]]==0) { ipuf--; }
1078         dx=node[arrayIndx[i]].nx-node[arrayIndx[near_node[ipuf]]].nx;
1079         dy=node[arrayIndx[i]].ny-node[arrayIndx[near_node[ipuf]]].ny;
1080         dz=node[arrayIndx[i]].nz-node[arrayIndx[near_node[ipuf]]].nz;
1081         if((abs(dx)<=local_gtol) &&(abs(dy)<=local_gtol) &&(abs(dz)<=local_gtol)) { n_closest_nodes*=2; goto repeatNodSearch; }
1082       }
1083       n=near_node[0]=ipuf;
1084       for(j=ipuf; j>0; j--)
1085       {
1086         if(arrayIndx[near_node[j]]==0) { n--; near_node[j]=-1; continue; }
1087         dx=node[arrayIndx[i]].nx-node[arrayIndx[near_node[j]]].nx;
1088         dy=node[arrayIndx[i]].ny-node[arrayIndx[near_node[j]]].ny;
1089         dz=node[arrayIndx[i]].nz-node[arrayIndx[near_node[j]]].nz;
1090         if(near_node[j]==i) { n--; near_node[j]=-1; }
1091         else if((abs(dx)>local_gtol) ||(abs(dy)>local_gtol) ||(abs(dz)>local_gtol)) { n--; near_node[j]=-1; }
1092       }
1093       if(n)
1094       {
1095         if(printFlag) printf("merge %d from:%d\n", i+1, k);
1096         j = mergeNode( arrayIndx, k, i, lock, near_node, nod2elem);
1097       }
1098     }
1099     if(rsort) free(rsort);
1100     if(nod2elem) { for(i=0; i<k; i++) free(nod2elem[i]); free(nod2elem); }
1101     if(arrayIndx) free(arrayIndx);
1102     if(orig_x) free( orig_x );
1103     if(orig_y) free( orig_y );
1104     if(orig_z) free( orig_z );
1105     if(sort_x) free( sort_x );
1106     if(sort_y) free( sort_y );
1107     if(sort_z) free( sort_z );
1108     if(sort_nx) free(sort_nx);
1109     if(sort_ny) free(sort_ny);
1110     if(sort_nz) free(sort_nz);
1111 
1112     makeSurfaces();
1113 
1114     /* when node coordinates were changed to the deformed ones then switch back  */
1115     if(addDispFlagLocal==2)
1116     {
1117       addDispToCoordinates(node);
1118     }
1119     else
1120     {
1121       realloc_colNr();
1122       updateDispLists();
1123     }
1124   }
1125   else if (( type[0]=='p' )||( type[0]=='P' ))
1126   {
1127     /* calculate all absolute r of all points and sort the indexes according to r */
1128     if ( (rsort = (Rsort *)malloc( (set[setNr].anz_p+1) * sizeof(Rsort))) == NULL )
1129       printf("ERROR: realloc failed: Rsort\n\n" );
1130     k=0;
1131     for( i=0; i<set[setNr].anz_p; i++)
1132     {
1133       rsort[k].r=sqrt(point[set[setNr].pnt[i]].px*point[set[setNr].pnt[i]].px+point[set[setNr].pnt[i]].py*point[set[setNr].pnt[i]].py+point[set[setNr].pnt[i]].pz*point[set[setNr].pnt[i]].pz);
1134       rsort[k].i=set[setNr].pnt[i];
1135       k++;
1136     }
1137     qsort( rsort, k, sizeof(Rsort), (void *)compareRsort );
1138 #if TEST
1139     for (i=0; i<k; i++)
1140       printf("%d p:%d point:%s r:%lf\n", i, rsort[i].i, point[rsort[i].i].name, rsort[i].r);
1141 #endif
1142 
1143     for(i=0; i<k; i++)
1144     {
1145       j = mergePnt( rsort,k, i, local_gtol, lock );
1146     }
1147     free(rsort);
1148 
1149     /* recalculate the line-shapes */
1150     for (i=0; i<anzGeo->l; i++) repLine(i);
1151     /* for (i=0; i<anzGeo->s; i++) repSurf(i);  */
1152   }
1153   else if (( type[0]=='e' )||( type[0]=='E' ))
1154   {
1155     if(addDispFlag==1)
1156     {
1157       addDispToCoordinates(node);
1158       // remember to switch back
1159       addDispFlagLocal=2;
1160     }
1161     mergeElem(setNr);
1162     makeSurfaces();
1163     realloc_colNr();
1164     if(addDispFlagLocal==2)
1165     {
1166       addDispToCoordinates(node);
1167     }
1168     else
1169     {
1170       updateDispLists();
1171     }
1172   }
1173   else if (( type[0]=='l' )||( type[0]=='L' ))
1174   {
1175     if((indx=(int *)malloc((set[setNr].anz_l+1)*sizeof(int)))==NULL)
1176     { errMsg("\nERROR: malloc failure in merge\n");
1177     return; }
1178     for (i=0; i<set[setNr].anz_l; i++)
1179     {
1180       indx[i]=set[setNr].line[i];
1181     }
1182     sum=i;
1183     for (i=0; i<sum; i++)
1184     {
1185       if (indx[i]==-1)
1186       {
1187         errMsg("ERROR: line is undefined, check program-code\n" );
1188       }
1189       else
1190       {
1191         mergeLine( indx[i], lock );
1192       }
1193     }
1194   }
1195   else if (( type[0]=='c' )||( type[0]=='C' ))
1196   {
1197     if((indx=(int *)malloc((set[setNr].anz_c+1)*sizeof(int)))==NULL)
1198     { errMsg("\nERROR: malloc failure in merge\n");
1199     return; }
1200     for (i=0; i<set[setNr].anz_c; i++)
1201     {
1202       indx[i]=set[setNr].lcmb[i];
1203     }
1204     sum=i;
1205     for (i=0; i<sum; i++)
1206     {
1207       if (indx[i]==-1)
1208       {
1209         errMsg("ERROR: lcmb is undefined, check program-code\n" );
1210       }
1211       else
1212       {
1213         mergeLcmb( indx[i] );
1214       }
1215     }
1216   }
1217   else if (( type[0]=='s' )||( type[0]=='S' ))
1218   {
1219     /* cyrcle through all surfs and add all lcmbs */
1220     for (i=0; i<set[setNr].anz_s; i++)
1221     {
1222       s= set[setNr].surf[i];
1223       for (j=0; j<surf[s].nl; j++)
1224       {
1225         l=surf[s].l[j];
1226         if (surf[s].typ[j]=='c')
1227         {
1228           seta( setNr, "c", l );
1229         }
1230       }
1231     }
1232     if((indx=(int *)malloc((set[setNr].anz_c+1)*sizeof(int)))==NULL)
1233     { errMsg("\nERROR: malloc failure in merge\n");
1234     return; }
1235     for (i=0; i<set[setNr].anz_c; i++)
1236     {
1237       indx[i]=set[setNr].lcmb[i];
1238     }
1239     sum=i;
1240     for (i=0; i<sum; i++)
1241     {
1242       if (indx[i]==-1)
1243       {
1244         errMsg("ERROR: lcmb is undefined, check program-code\n" );
1245       }
1246       else
1247       {
1248         j = mergeLcmb( indx[i] );
1249       }
1250     }
1251     free(indx);
1252     if((indx=(int *)malloc((set[setNr].anz_s+1)*sizeof(int)))==NULL)
1253     { errMsg("\nERROR: malloc failure in merge\n");
1254     return; }
1255     for (i=0; i<set[setNr].anz_s; i++)
1256     {
1257       indx[i]=set[setNr].surf[i];
1258     }
1259     sum=i;
1260     for (i=0; i<sum; i++)
1261     {
1262       if (indx[i]==-1)
1263       {
1264         errMsg("ERROR: surf is undefined, check program-code\n" );
1265       }
1266       else
1267       {
1268         j = mergeSurf( indx[i], local_gtol );
1269       }
1270     }
1271   }
1272   else
1273   {
1274     printf (" entity:%c not known\n", type[0]);
1275     return;
1276   }
1277   if (indx!=NULL) free(indx); indx=NULL;
1278   printf("ready\n");
1279 }
1280 
1281 
1282