FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hecmw_dlb_elem_sr_adapt.c
Go to the documentation of this file.
1 /*****************************************************************************
2  * Copyright (c) 2019 FrontISTR Commons
3  * This software is released under the MIT License, see LICENSE.txt
4  *****************************************************************************/
5 
6 #include "hecmw_repart.h"
7 
8 extern struct hecmwST_local_mesh *mesh;
9 extern struct hecmwST_result_data *data;
10 extern struct hecmwST_local_mesh *new_mesh;
11 extern struct hecmwST_result_data *new_data;
12 extern int stack_part_send_recv(int neibpetot, int *neibpe, int *stack_import,
13  int *stack_export, HECMW_Comm repart_comm,
14  int my_rank);
15 extern int stack_whole_send_recv(int pesize, int *stack_import,
16  int *stack_export, HECMW_Comm repart_comm,
17  int my_rank);
18 
19 extern int int_part_send_recv(int n, int neibpetot, int *neibpe,
20  int *stack_import, int *nod_import,
21  int *stack_export, int *nod_export, int *x,
22  HECMW_Comm repart_comm, int my_rank);
23 extern int double_part_send_recv(int n, int neibpetot, int *neibpe,
24  int *stack_import, int *nod_import,
25  int *stack_export, int *nod_export, double *x,
26  HECMW_Comm repart_comm, int my_rank);
27 extern void int_whole_send_recv(int n1, int n2, int pesize, int *stack_import,
28  int *nod_import, int *stack_export,
29  int *nod_export, int *x, int *y,
30  HECMW_Comm repart_comm, int my_rank);
31 extern void double_whole_send_recv(int n1, int n2, int pesize,
32  int *stack_import, int *nod_import,
33  int *stack_export, int *nod_export,
34  double *x, double *y, HECMW_Comm repart_comm,
35  int my_rank);
36 
37 extern void int2_whole_send_recv(int n1, int n2, int pesize, int *stack_import,
38  int *stack_export, int *x, int *y,
39  HECMW_Comm repart_comm, int my_rank);
40 extern void int3_whole_send_recv(int n1, int n2, int pesize, int *stack_import,
41  int *stack_export, int *x, int *y,
42  HECMW_Comm repart_comm, int my_rank);
43 extern void double2_whole_send_recv(int n1, int n2, int pesize,
44  int *stack_import, int *stack_export,
45  double *x, double *y,
46  HECMW_Comm repart_comm, int my_rank);
47 extern void whole_copy_array(int *recv_elem_num, int *global_recv_elem_num,
48  int mynode, int pesize, HECMW_Comm repart_comm);
49 
50 void mesh_migration_adapt(int mynode, int pesize, Result_part *result,
51  int *vtxdist)
52 
53 {
54  int *send_elem_num, *send_elem, *recv_elem_num, *recv_elem, *recv_parent_num,
55  *count_elem, *send_parent, *send_parent_num, *count_parent;
56  int *flag_hit, *flag_elem_hit;
57  int ncon, id_elem;
58  int i, j, k, m;
59  int *tmp_int_send, *tmp_int_recv, *tmp2_int_send, *tmp2_int_recv;
60  int *new_part;
61  int *tmp_stack_send, *tmp_stack_recv, *tmp_neibpe, tmp_int, tmp_sum;
62  int *send_ptr_num, *send_ptr_parent_num, *recv_ptr_num, *recv_ptr_parent_num;
63  int *global_index, *global_index_hit;
64  int *send_node_num, *recv_node_num, *send_node, *recv_node;
65  int *tmp_node, *count_node, *tmp_send, *tmp_recv;
66  int *import_index, *export_index, *count_index, *export_node, *import_node;
67  int import_n_neighbor_pe, export_n_neighbor_pe, *import_neighbor_pe,
68  *export_neighbor_pe;
69  /* int *import2_index, *export2_index, *export2_node, *import2_node;
70  */
71  int *tmp_int_nodeid, *tmp2_int_nodeid;
72  double *tmp_node_d, *recv_node_d, *tmp2_node_d, *recv2_node_d, *tmp_send_d,
73  *tmp_recv_d, *tmp2_send_d, *tmp2_recv_d;
74  int *tmp_node_i, *recv_node_i, *tmp2_node_i, *recv2_node_i;
75  HECMW_Status stat;
76  FILE *test_fp, *test_fp2;
77  char test_file[128];
78  int min_pe, tmp_pe, local_nid;
79  int nvtxs, tmp_nvtxs, tmp_peid, tmp_lid;
80  int *new_elem, *global_recv_elem_num, *global_recv_parent_num;
81  int *send_adapt_parent, *recv_adapt_parent, *send_tmp, *recv_tmp,
82  *send_tmp_num, *recv_tmp_num, *send_adapt_ptype, *count_num;
83  int t_elem, l_child, new_l_child;
84  int *inter_elem, *send_adapt_child, *recv_adapt_child, *send_adapt_child_num,
85  *recv_adapt_child_num, *send_index_child, *recv_index_child;
86  int *send2_adapt_child_num, *recv2_adapt_child_num;
87  int *send_inter, *recv_inter, *send_inter_num, *recv_inter_num,
88  *global_recv_inter_num;
89  int *send_parent_inter, *recv_parent_inter, *send_parent_inter_num,
90  *recv_parent_inter_num, *global_recv_parent_inter_num;
91  int init_flag, new_elemid;
92  int *tmp_elem_grp, num_grp_item, *tmp_surf_grp, *tmp_surf_id, *tmp2_surf_id;
93  Tmp_grp_inf *tmp_grp;
94  int *new_nelem_dist;
95  int tn_component;
96  int *new_vtxdist;
97  int *new2old, *count_elem_index, *new_tmp, *new_tmp2, *old2new;
98  double t1, t2;
99 
100  if (mynode == 0) {
101  t1 = HECMW_Wtime();
102  }
103 
104 #ifdef TEST
105  sprintf(test_file, "test3.%d", mynode);
106  test_fp = fopen(test_file, "w");
107  if (test_fp == NULL) {
108  fprintf(stderr, "Cannot open test_file\n");
109  exit(0);
110  }
111 #endif
112 /* fprintf(test_fp, "%d\n", mesh->ne_internal);
113  for(i=0;i<mesh->ne_internal;i++)
114  fprintf(test_fp, "%d %d\n", mesh->elem_internal_list[i],
115  mesh->adapt_type[mesh->elem_internal_list[i]]);
116  fprintf(test_fp, "%d %d\n", mesh->n_node, mesh->nn_internal);
117 
118  fprintf(test_fp, "%d \n", mesh->n_neighbor_pe);
119  for(i=0;i<mesh->n_neighbor_pe;i++)
120  fprintf(test_fp, "%d ", mesh->neighbor_pe[i]);
121  fprintf(test_fp, "\n");
122  for(i=0;i<mesh->n_neighbor_pe+1;i++)
123  fprintf(test_fp, "%d ", mesh->import_index[i]);
124  fprintf(test_fp, "\n");
125  for(i=0;i<mesh->n_neighbor_pe+1;i++)
126  fprintf(test_fp, "%d ", mesh->export_index[i]);
127  fprintf(test_fp, "\n");
128  for(i=0;i<mesh->import_index[mesh->n_neighbor_pe];i++)
129  fprintf(test_fp, "%d\n", mesh->import_item[i]);
130  fprintf(test_fp, "export node*****\n");
131  for(i=0;i<mesh->export_index[mesh->n_neighbor_pe];i++)
132  fprintf(test_fp, "%d\n", mesh->export_item[i]);
133 */
134 /*
135  fprintf(test_fp, "in PE %d n_elem=%d ne_internal=%d\n", mynode,
136  mesh->n_elem, mesh->ne_internal);
137  for(i=0;i<mesh->n_elem;i++)
138  fprintf(test_fp, "i=%d %d %d %d %d %d %d\n", i, mesh->elem_ID[i*2],
139  mesh->elem_ID[i*2+1], mesh->adapt_parent[i*2], mesh->adapt_parent[i*2+1],
140  mesh->adapt_type[i], mesh->elem_node_index[i]);
141  fclose(test_fp);
142 */
143 #ifdef test
144  sprintf(test_file, "test5.%d", mynode);
145  test_fp2 = fopen(test_file, "w");
146  if (test_fp == NULL) HECMW_dlb_print_exit("Cannot open test_file\n");
147 #endif
148 
149  /* 1: first update new partition to include import nodes */
150 
151  new_part = (int *)calloc(mesh->n_node, sizeof(int));
152  if (new_part == NULL) HECMW_dlb_memory_exit("new_part");
153  for (i = 0; i < mesh->nn_internal; i++) new_part[i] = result->part[i];
156  mesh->export_item, new_part, mesh->HECMW_COMM, mynode);
157  free(result->part);
158  result->part = new_part;
159  new_elem = (int *)calloc(mesh->n_elem, sizeof(int));
160  if (new_elem == NULL) HECMW_dlb_memory_exit("new_elem");
161  inter_elem = (int *)calloc(mesh->n_elem, sizeof(int));
162  /* tmp_int=0;
163  */
164  for (i = 0; i < mesh->n_elem; i++) {
165  inter_elem[i] = -1;
166  if (mesh->elem_ID[i * 2 + 1] == mynode) {
167  inter_elem[mesh->elem_ID[i * 2] - 1] = i;
168  /* tmp_int++;
169  */
170  }
171  }
172  /* fprintf(test_fp, "ne_internal====%d\n", mesh->ne_internal);
173  for(i=0;i<mesh->n_elem;i++)
174  fprintf(test_fp, "%d %d %d %d\n", i, inter_elem[i],
175  mesh->elem_ID[i*2], mesh->elem_ID[i*2+1]);
176  fclose(test_fp);
177  */
178  for (i = 0; i < mesh->n_elem; i++) new_elem[i] = -1;
179  if (pesize > 1)
181  mesh->HECMW_COMM);
182  else
183  t_elem = mesh->n_elem;
184 
185  /* send_parent_num=(int *)calloc(pesize+1, sizeof(int)); */
186  send_elem_num = (int *)calloc(pesize + 1, sizeof(int));
187  send_inter_num = (int *)calloc(pesize + 1, sizeof(int));
188  recv_elem_num = (int *)calloc(pesize + 1, sizeof(int));
189  recv_inter_num = (int *)calloc(pesize + 1, sizeof(int));
190  if ((send_elem_num == NULL) || (send_inter_num == NULL) ||
191  (recv_elem_num == NULL) || (recv_inter_num == NULL))
192  HECMW_dlb_memory_exit("send_elem_num ");
193  flag_hit = (int *)calloc(pesize, sizeof(int));
194  flag_elem_hit = (int *)calloc(mesh->n_elem, sizeof(int));
195  if ((flag_hit == NULL) || (flag_elem_hit == NULL))
196  HECMW_dlb_memory_exit("flag_hit");
197 
198  for (i = 0; i < pesize + 1; i++) {
199  send_elem_num[i] = 0;
200  send_inter_num[i] = 0;
201  recv_elem_num[i] = 0;
202  recv_inter_num[i] = 0;
203  }
204 
205  for (i = 0; i < mesh->n_elem; i++) {
206  for (j = 0; j < pesize; j++) flag_hit[j] = 0;
207  init_flag = pesize;
208  if ((mesh->elem_ID[i * 2 + 1] == mynode) && (mesh->adapt_type[i] == 0)) {
209  ncon = mesh->elem_node_index[i + 1] - mesh->elem_node_index[i];
210  for (j = mesh->elem_node_index[i]; j < mesh->elem_node_index[i + 1];
211  j++) {
212  if (flag_hit[result->part[mesh->elem_node_item[j] - 1]] == 0) {
213  flag_hit[result->part[mesh->elem_node_item[j] - 1]] = 1;
214  send_elem_num[result->part[mesh->elem_node_item[j] - 1] + 1]++;
215  if (result->part[mesh->elem_node_item[j] - 1] < init_flag)
216  init_flag = result->part[mesh->elem_node_item[j] - 1];
217  }
218  }
219  send_inter_num[init_flag + 1]++;
220  }
221  }
222  /* for(i=1;i<pesize+1;i++)
223  fprintf(stderr, "in PE %d the number send to PE %d is
224  %d\n", mynode, i, send_elem_num[i]);
225  */
226 
227  for (i = 1; i < pesize + 1; i++) {
228  send_elem_num[i] = send_elem_num[i - 1] + send_elem_num[i];
229  send_inter_num[i] = send_inter_num[i - 1] + send_inter_num[i];
230  }
231  /*
232  for(i=0;i<pesize+1;i++)
233  fprintf(stderr, "PE %d: send %d is %d\n", mynode, i,
234  send_elem_num[i]);
235  */
237 
238  stack_whole_send_recv(pesize, send_elem_num, recv_elem_num, mesh->HECMW_COMM,
239  mynode);
240  stack_whole_send_recv(pesize, send_inter_num, recv_inter_num,
241  mesh->HECMW_COMM, mynode);
242  /*
243  for(i=0;i<pesize+1;i++)
244  fprintf(stderr, "PE %d: recv %d is %d\n", mynode, i,
245  recv_elem_num[i]);
246  */
247 
248  global_recv_elem_num = (int *)calloc(pesize * (pesize + 1), sizeof(int));
249  global_recv_inter_num = (int *)calloc(pesize * (pesize + 1), sizeof(int));
250  if ((global_recv_elem_num == NULL) || (global_recv_inter_num == NULL))
251  HECMW_dlb_memory_exit("global_recv_elem_num");
252  whole_copy_array(recv_elem_num, global_recv_elem_num, mynode, pesize,
253  mesh->HECMW_COMM);
254  whole_copy_array(recv_inter_num, global_recv_inter_num, mynode, pesize,
255  mesh->HECMW_COMM);
256  /* for(i=1;i<pesize+1;i++)
257  fprintf(stderr, "AFTER--in PE0 the number send to PE
258  %d is %d\n", i, send_elem_num[i]);
259  */
260  if (send_elem_num[pesize] > 0) {
261  send_elem = (int *)calloc(send_elem_num[pesize], sizeof(int));
262  count_elem = (int *)calloc(pesize, sizeof(int));
263  count_num = (int *)calloc(pesize, sizeof(int));
264  send_inter = (int *)calloc(send_elem_num[pesize], sizeof(int));
265 
266  if ((send_elem == NULL) || (count_elem == NULL))
267  HECMW_dlb_memory_exit("send_elem and count_elem");
268  for (i = 0; i < pesize; i++) {
269  count_elem[i] = 0;
270  count_num[i] = 0;
271  }
272  }
273  for (i = 0; i < mesh->n_elem; i++) {
274  for (j = 0; j < pesize; j++) flag_hit[j] = 0;
275  init_flag = pesize;
276  if ((mesh->elem_ID[i * 2 + 1] == mynode) && (mesh->adapt_type[i] == 0)) {
277  ncon = mesh->elem_node_index[i + 1] - mesh->elem_node_index[i];
278  for (j = mesh->elem_node_index[i]; j < mesh->elem_node_index[i + 1];
279  j++) {
280  if (flag_hit[result->part[mesh->elem_node_item[j] - 1]] == 0) {
281  flag_hit[result->part[mesh->elem_node_item[j] - 1]] = 1;
282  send_elem[send_elem_num[result->part[mesh->elem_node_item[j] - 1]] +
283  count_elem[result->part[mesh->elem_node_item[j] - 1]]] = i;
284  count_elem[result->part[mesh->elem_node_item[j] - 1]]++;
285  if (result->part[mesh->elem_node_item[j] - 1] < init_flag)
286  init_flag = result->part[mesh->elem_node_item[j] - 1];
287  }
288  }
289  new_elemid = init_flag * t_elem +
290  global_recv_inter_num[init_flag * (pesize + 1) + mynode] +
291  count_num[init_flag];
292  new_elem[i] = new_elemid;
293  /* new_elemid=init_flag*t_elem+global_recv_elem_num[init_flag*(pesize+1)+mynode]+count_elem[init_flag]-1;
294  new_elem[i]=new_elemid;
295  */
296  count_num[init_flag]++;
297  }
298  }
299  for (i = 0; i < send_elem_num[pesize]; i++)
300  send_inter[i] = new_elem[send_elem[i]];
301 
302  /* -------start finding parent information ----------------*/
303  /* find all the parents which need be sent */
304  send_parent_num = (int *)calloc(pesize + 1, sizeof(int));
305  send_parent_inter_num = (int *)calloc(pesize + 1, sizeof(int));
306  recv_parent_num = (int *)calloc(pesize + 1, sizeof(int));
307  recv_parent_inter_num = (int *)calloc(pesize + 1, sizeof(int));
308 
309  for (i = 0; i < pesize + 1; i++) {
310  send_parent_num[i] = 0;
311  send_parent_inter_num[i] = 0;
312  }
313  for (i = 0; i < mesh->n_elem; i++) flag_elem_hit[i] = 0;
314  for (j = 0; j < pesize; j++) {
315  send_parent_num[j + 1] = 0;
316  for (k = send_elem_num[j]; k < send_elem_num[j + 1]; k++) {
317  id_elem = send_elem[k];
318  while ((mesh->adapt_parent[id_elem * 2] > 0) &&
319  (mesh->adapt_parent[id_elem * 2 + 1] == mynode) &&
320  (flag_elem_hit[mesh->adapt_parent[id_elem * 2] - 1] == 0)) {
321  flag_elem_hit[mesh->adapt_parent[id_elem * 2] - 1] = 1;
322  send_parent_num[j + 1]++;
323  id_elem = inter_elem[mesh->adapt_parent[id_elem * 2] - 1];
324  }
325  }
326  }
327 
328  for (i = 1; i < pesize + 1; i++) {
329  send_parent_num[i] = send_parent_num[i - 1] + send_parent_num[i];
330  }
331 
332  recv_parent_num = (int *)calloc(pesize + 1, sizeof(int));
333  if (recv_parent_num == NULL) HECMW_dlb_memory_exit("recv_parent_num");
334  stack_whole_send_recv(pesize, send_parent_num, recv_parent_num,
335  mesh->HECMW_COMM, mynode);
336  global_recv_parent_num = (int *)calloc(pesize * (pesize + 1), sizeof(int));
337  if (global_recv_parent_num == NULL)
338  HECMW_dlb_memory_exit("global_recv_parent_num");
339  whole_copy_array(recv_parent_num, global_recv_parent_num, mynode, pesize,
340  mesh->HECMW_COMM);
341 
342  send_parent = (int *)calloc(send_parent_num[pesize], sizeof(int));
343  send_parent_inter = (int *)calloc(send_parent_num[pesize], sizeof(int));
344  count_parent = (int *)calloc(pesize, sizeof(int));
345 
346  if ((send_parent == NULL) || (count_parent == NULL))
347  HECMW_dlb_memory_exit("send_parent and count_parent");
348  /* to keep each parent just occur once in all PEs */
349  for (i = 0; i < mesh->n_elem; i++) flag_elem_hit[i] = 0;
350  for (j = 0; j < pesize; j++) {
351  count_parent[j] = 0;
352  for (k = send_elem_num[j]; k < send_elem_num[j + 1]; k++) {
353  id_elem = send_elem[k];
354  while ((mesh->adapt_parent[id_elem * 2] > 0) &&
355  (mesh->adapt_parent[id_elem * 2 + 1] == mynode) &&
356  (flag_elem_hit[mesh->adapt_parent[id_elem * 2] - 1] == 0)) {
357  flag_elem_hit[mesh->adapt_parent[id_elem * 2] - 1] = 1;
358  send_parent[send_parent_num[j] + count_parent[j]] =
359  inter_elem[mesh->adapt_parent[id_elem * 2] - 1];
360  /* send_parent_inter[send_parent_num[j]+count_parent[j]]=j*t_elem+global_recv_inter_num[j*(pesize+1)+pesize]+global_recv_parent_num[j*(pesize+1)+
361  mynode]+count_parent[j];
362  */
363 
364  count_parent[j]++;
365 
366  new_elem[inter_elem[mesh->adapt_parent[id_elem * 2] - 1]] =
367  j * t_elem + global_recv_inter_num[j * (pesize + 1) + pesize] +
368  global_recv_parent_num[j * (pesize + 1) + mynode] +
369  count_parent[j] - 1;
370 
371  id_elem = inter_elem[mesh->adapt_parent[id_elem * 2] - 1];
372  }
373  }
374  }
375  for (i = 0; i < send_parent_num[pesize]; i++) {
376  send_parent_inter[i] = new_elem[send_parent[i]];
377  if (send_parent_inter[i] < 0)
378  HECMW_dlb_print_exit("There is something wrong with send_parent_inter");
379  }
380  /* ------------------- start sending parent information ----------------*/
381 
382  send_tmp_num = (int *)calloc(pesize + 1, sizeof(int));
383  recv_tmp_num = (int *)calloc(pesize + 1, sizeof(int));
384  if ((send_tmp_num == NULL) || (recv_tmp_num == NULL))
385  HECMW_dlb_memory_exit("send_tmp_num");
386  if (send_elem_num[pesize] > 0) {
387  send_adapt_parent = (int *)calloc(send_elem_num[pesize], sizeof(int));
388  send_adapt_ptype = (int *)calloc(send_elem_num[pesize], sizeof(int));
389  if ((send_adapt_parent == NULL) || (send_adapt_ptype == NULL))
390  HECMW_dlb_memory_exit("send_adapt_parent");
391  }
392 
393  for (i = 0; i < pesize + 1; i++) send_tmp_num[i] = 0;
394 
395  for (i = 0; i < send_elem_num[pesize]; i++) {
396  if (mesh->adapt_parent[send_elem[i] * 2 + 1] < 0)
397  send_adapt_parent[i] = -1;
398  else if (mesh->adapt_parent[send_elem[i] * 2 + 1] != mynode)
399  send_tmp_num[mesh->adapt_parent[send_elem[i] * 2 + 1] + 1]++;
400  else if (new_elem[inter_elem[mesh->adapt_parent[send_elem[i] * 2] - 1]] !=
401  -1)
402  send_adapt_parent[i] =
403  new_elem[inter_elem[mesh->adapt_parent[send_elem[i] * 2] - 1]];
404  else if (new_elem[inter_elem[mesh->adapt_parent[send_elem[i] * 2] - 1]] ==
405  -1) {
406  fprintf(stderr, "There is something wrong with parent information\n");
407  fprintf(stderr, "i=%d, send_elem[i]=%d parent is %d PE=%d\n", i,
408  send_elem[i], mesh->adapt_parent[send_elem[i] * 2] - 1,
409  mesh->adapt_parent[send_elem[i] * 2]);
410  HECMW_dlb_print_exit("Error in parent finding");
411  }
412  }
413  for (i = 1; i < pesize + 1; i++) {
414  send_tmp_num[i] = send_tmp_num[i - 1] + send_tmp_num[i];
415  }
416 
417  stack_whole_send_recv(pesize, send_tmp_num, recv_tmp_num, mesh->HECMW_COMM,
418  mynode);
419 
420  for (i = 0; i < pesize; i++) count_num[i] = 0;
421  if (send_tmp_num[pesize] > 0) {
422  send_tmp = (int *)calloc(send_tmp_num[pesize], sizeof(int));
423  if (send_tmp == NULL) HECMW_dlb_memory_exit("send_tmp");
424  }
425  if (recv_tmp_num[pesize] > 0) {
426  recv_tmp = (int *)calloc(recv_tmp_num[pesize], sizeof(int));
427  if (recv_tmp == NULL) HECMW_dlb_memory_exit("send_tmp, recv_tmp");
428  }
429  for (i = 0; i < send_elem_num[pesize]; i++) {
430  if (mesh->adapt_parent[send_elem[i] * 2 + 1] < 0)
431  send_adapt_parent[i] = -1;
432  else if (mesh->adapt_parent[send_elem[i] * 2 + 1] != mynode) {
433  send_tmp[send_tmp_num[mesh->adapt_parent[send_elem[i] * 2 + 1]] +
434  count_num[mesh->adapt_parent[send_elem[i] * 2 + 1]]] =
435  inter_elem[mesh->adapt_parent[send_elem[i] * 2] - 1];
436  count_num[mesh->adapt_parent[send_elem[i] * 2 + 1]]++;
437  }
438  }
439 
440  if (send_tmp_num[pesize] > 0)
441  int2_whole_send_recv(send_tmp_num[pesize], recv_tmp_num[pesize], pesize,
442  recv_tmp_num, send_tmp_num, send_tmp, recv_tmp,
443  mesh->HECMW_COMM, mynode);
445  for (i = 0; i < recv_tmp_num[pesize]; i++) {
446  if (new_elem[recv_tmp[i]] == -1)
447  HECMW_dlb_print_exit("There is something wrong in parent inf");
448  else
449  recv_tmp[i] = new_elem[recv_tmp[i]];
450  }
451  if (send_tmp_num[pesize] > 0)
452  int2_whole_send_recv(recv_tmp_num[pesize], send_tmp_num[pesize], pesize,
453  send_tmp_num, recv_tmp_num, recv_tmp, send_tmp,
454  mesh->HECMW_COMM, mynode);
455 
456  for (i = 0; i < pesize; i++) count_num[i] = 0;
457  for (i = 0; i < send_elem_num[pesize]; i++) {
458  if (mesh->adapt_parent[send_elem[i] * 2 + 1] < 0) {
459  send_adapt_parent[i] = -1;
460  send_adapt_ptype[i] = -1;
461  } else if (mesh->adapt_parent[send_elem[i] * 2 + 1] != mynode) {
462  send_adapt_parent[i] =
463  send_tmp[send_tmp_num[mesh->adapt_parent[send_elem[i] * 2 + 1]] +
464  count_num[mesh->adapt_parent[send_elem[i] * 2 + 1]]];
465  send_adapt_ptype[i] = mesh->adapt_parent_type[send_elem[i]];
466 
467  count_num[mesh->adapt_parent[send_elem[i] * 2 + 1]]++;
468  } else if (new_elem[inter_elem[mesh->adapt_parent[send_elem[i] * 2] - 1]] !=
469  -1) {
470  send_adapt_parent[i] =
471  new_elem[inter_elem[mesh->adapt_parent[send_elem[i] * 2] - 1]];
472  send_adapt_ptype[i] = mesh->adapt_parent_type[send_elem[i]];
473  } else if (new_elem[inter_elem[mesh->adapt_parent[send_elem[i] * 2] - 1]] ==
474  -1)
475  HECMW_dlb_print_exit("There is something wrong with parent information");
476  }
477  if (recv_elem_num[pesize] > 0) {
478  tmp_int_recv = (int *)calloc(recv_elem_num[pesize], sizeof(int));
479  if (tmp_int_recv == NULL) HECMW_dlb_memory_exit("tmp_int_recv");
480  int2_whole_send_recv(send_elem_num[pesize], recv_elem_num[pesize], pesize,
481  recv_elem_num, send_elem_num, send_adapt_parent,
482  tmp_int_recv, mesh->HECMW_COMM, mynode);
483  }
484  /* new mesh whole inf. copy */
485  /* for(i=0;i<mesh->n_elem;i++) {
486  if(new_elem[i]==-1)
487  fprintf(stderr, "Wrong in new elem: %d %d\n", mynode, i);
488  HECMW_dlb_print_exit("check");
489  }
490  */
492  new_mesh->files = (char **)calloc(mesh->hecmw_n_file, sizeof(char *));
493  for (i = 9; i < mesh->hecmw_n_file; i++)
494  new_mesh->files[i] = (char *)calloc(128, sizeof(char));
495 
496  sprintf(new_mesh->header, "%s", mesh->header);
497  sprintf(new_mesh->gridfile, "%s", mesh->gridfile);
504  new_mesh->n_dof = mesh->n_dof;
506 
507  new_mesh->n_elem = recv_elem_num[pesize] + recv_parent_num[pesize];
510  if (new_mesh->n_elem <= 0) HECMW_dlb_print_exit("Error: New mesh: n_elem==0");
511  new_mesh->adapt_parent = (int *)calloc(new_mesh->n_elem * 2, sizeof(int));
512  if (new_mesh->adapt_parent == NULL)
513  HECMW_dlb_memory_exit("new_mesh: adaption_parent");
514  new_mesh->adapt_parent_type = (int *)calloc(new_mesh->n_elem, sizeof(int));
516  HECMW_dlb_memory_exit("new_mesh: adapt_parent_type");
517  for (i = 0; i < recv_elem_num[pesize]; i++) {
518  if (tmp_int_recv[i] == -1) {
519  new_mesh->adapt_parent[i * 2] = 0;
520  new_mesh->adapt_parent[i * 2 + 1] = -1;
521  } else {
522  new_mesh->adapt_parent[i * 2] = (tmp_int_recv[i] % t_elem) + 1;
523  new_mesh->adapt_parent[i * 2 + 1] = tmp_int_recv[i] / t_elem;
524  }
525  }
526  int2_whole_send_recv(send_elem_num[pesize], recv_elem_num[pesize], pesize,
527  recv_elem_num, send_elem_num, send_adapt_ptype,
528  tmp_int_recv, mesh->HECMW_COMM, mynode);
529  for (i = 0; i < recv_elem_num[pesize]; i++)
530  new_mesh->adapt_parent_type[i] = tmp_int_recv[i];
531  if (send_tmp_num[pesize] > 0) free(send_tmp);
532  if (recv_tmp_num[pesize] > 0) free(recv_tmp);
533  if (send_elem_num[pesize] > 0) {
534  free(send_adapt_parent);
535  free(send_adapt_ptype);
536  }
537  if (send_parent_num[pesize] > 0) {
538  send_adapt_parent = (int *)calloc(send_parent_num[pesize], sizeof(int));
539  send_adapt_ptype = (int *)calloc(send_parent_num[pesize], sizeof(int));
540  if ((send_adapt_parent == NULL) || (send_adapt_ptype == NULL))
541  HECMW_dlb_memory_exit("send_adapt_parent");
542  }
543 
544  for (i = 0; i < pesize + 1; i++) send_tmp_num[i] = 0;
545 
546  for (i = 0; i < send_parent_num[pesize]; i++) {
547  if (mesh->adapt_parent[send_parent[i] * 2 + 1] < 0)
548  send_adapt_parent[i] = -1;
549  else if (mesh->adapt_parent[send_parent[i] * 2 + 1] != mynode)
550  send_tmp_num[mesh->adapt_parent[send_parent[i] * 2 + 1] + 1]++;
551  else if (new_elem[inter_elem[mesh->adapt_parent[send_parent[i] * 2] - 1]] !=
552  -1)
553  send_adapt_parent[i] =
554  new_elem[inter_elem[mesh->adapt_parent[send_parent[i] * 2] - 1]];
555  else if (new_elem[inter_elem[mesh->adapt_parent[send_parent[i] * 2] - 1]] ==
556  -1)
557  HECMW_dlb_print_exit("There is something wrong with parent information");
558  }
559  for (i = 1; i < pesize + 1; i++) {
560  send_tmp_num[i] = send_tmp_num[i - 1] + send_tmp_num[i];
561  }
562 
563  stack_whole_send_recv(pesize, send_tmp_num, recv_tmp_num, mesh->HECMW_COMM,
564  mynode);
565 
566  for (i = 0; i < pesize; i++) count_num[i] = 0;
567  if (send_tmp_num[pesize] > 0) {
568  send_tmp = (int *)calloc(send_tmp_num[pesize], sizeof(int));
569  if (send_tmp == NULL) HECMW_dlb_memory_exit("send_tmp, recv_tmp");
570  }
571  if (recv_tmp_num[pesize] > 0) {
572  recv_tmp = (int *)calloc(recv_tmp_num[pesize], sizeof(int));
573  if (recv_tmp == NULL) HECMW_dlb_memory_exit("send_tmp, recv_tmp");
574  }
575  for (i = 0; i < send_parent_num[pesize]; i++) {
576  if (mesh->adapt_parent[send_parent[i] * 2 + 1] < 0)
577  send_adapt_parent[i] = -1;
578  else if (mesh->adapt_parent[send_parent[i] * 2 + 1] != mynode) {
579  send_tmp[send_tmp_num[mesh->adapt_parent[send_parent[i] * 2 + 1]] +
580  count_num[mesh->adapt_parent[send_parent[i] * 2 + 1]]] =
581  inter_elem[mesh->adapt_parent[send_parent[i] * 2] - 1];
582  count_num[mesh->adapt_parent[send_parent[i] * 2 + 1]]++;
583  }
584  }
585  if (send_tmp_num[pesize] > 0)
586  int2_whole_send_recv(send_tmp_num[pesize], recv_tmp_num[pesize], pesize,
587  recv_tmp_num, send_tmp_num, send_tmp, recv_tmp,
588  mesh->HECMW_COMM, mynode);
589  for (i = 0; i < recv_tmp_num[pesize]; i++) {
590  if (new_elem[recv_tmp[i]] == -1)
591  HECMW_dlb_print_exit("There is something wrong in parent inf");
592  else
593  recv_tmp[i] = new_elem[recv_tmp[i]];
594  }
595  if (send_tmp_num[pesize] > 0)
596  int2_whole_send_recv(recv_tmp_num[pesize], send_tmp_num[pesize], pesize,
597  send_tmp_num, recv_tmp_num, recv_tmp, send_tmp,
598  mesh->HECMW_COMM, mynode);
599 
600  for (i = 0; i < pesize; i++) count_num[i] = 0;
601  for (i = 0; i < send_parent_num[pesize]; i++) {
602  if (mesh->adapt_parent[send_parent[i] * 2 + 1] < 0) {
603  send_adapt_parent[i] = -1;
604  send_adapt_ptype[i] = -1;
605  } else if (mesh->adapt_parent[send_parent[i] * 2 + 1] != mynode) {
606  send_adapt_parent[i] =
607  send_tmp[send_tmp_num[mesh->adapt_parent[send_parent[i] * 2 + 1]] +
608  count_num[mesh->adapt_parent[send_parent[i] * 2 + 1]]];
609  send_adapt_ptype[i] = mesh->adapt_parent_type[send_parent[i]];
610 
611  count_num[mesh->adapt_parent[send_parent[i] * 2 + 1]]++;
612  } else if (new_elem[inter_elem[mesh->adapt_parent[send_parent[i] * 2] -
613  1]] != -1) {
614  send_adapt_parent[i] =
615  new_elem[inter_elem[mesh->adapt_parent[send_parent[i] * 2] - 1]];
616  send_adapt_ptype[i] = mesh->adapt_parent_type[send_parent[i]];
617  } else if (new_elem[inter_elem[mesh->adapt_parent[send_parent[i] * 2] -
618  1]] == -1)
619  HECMW_dlb_print_exit("There is something wrong with parent information");
620  }
621 
622  tmp2_int_recv = (int *)calloc(recv_parent_num[pesize], sizeof(int));
623  if (tmp2_int_recv == NULL) HECMW_dlb_memory_exit("tmp2_int_recv");
624  int2_whole_send_recv(send_parent_num[pesize], recv_parent_num[pesize], pesize,
625  recv_parent_num, send_parent_num, send_adapt_parent,
626  tmp2_int_recv, mesh->HECMW_COMM, mynode);
627 
628  for (i = 0; i < recv_parent_num[pesize]; i++) {
629  if (tmp2_int_recv[i] == -1) {
630  new_mesh->adapt_parent[(i + recv_elem_num[pesize]) * 2] = 0;
631  new_mesh->adapt_parent[(i + recv_elem_num[pesize]) * 2 + 1] = -1;
632  } else {
633  new_mesh->adapt_parent[(i + recv_elem_num[pesize]) * 2] =
634  (tmp2_int_recv[i] % t_elem) + 1;
635  new_mesh->adapt_parent[(i + recv_elem_num[pesize]) * 2 + 1] =
636  tmp2_int_recv[i] / t_elem;
637  }
638  }
639  int2_whole_send_recv(send_parent_num[pesize], recv_parent_num[pesize], pesize,
640  recv_parent_num, send_parent_num, send_adapt_ptype,
641  tmp2_int_recv, mesh->HECMW_COMM, mynode);
642  for (i = 0; i < recv_parent_num[pesize]; i++)
643  new_mesh->adapt_parent_type[i + recv_elem_num[pesize]] = tmp2_int_recv[i];
644  if (send_tmp_num[pesize] > 0) free(send_tmp);
645  if (recv_tmp_num[pesize] > 0) free(recv_tmp);
646  if (send_parent_num[pesize] > 0) {
647  free(send_adapt_parent);
648  free(send_adapt_ptype);
649  }
651 
652  /* ------------- start migration child information ------------------ */
653  send_adapt_child_num = (int *)calloc(pesize + 1, sizeof(int));
654  for (j = 0; j < pesize + 1; j++) send_adapt_child_num[j] = 0;
655  for (j = 0; j < pesize; j++) {
656  for (i = send_elem_num[j]; i < send_elem_num[j + 1]; i++)
657  send_adapt_child_num[j + 1] +=
658  mesh->adapt_children_index[send_elem[i] + 1] -
659  mesh->adapt_children_index[send_elem[i]];
660  }
661  for (i = 1; i < pesize + 1; i++) {
662  send_adapt_child_num[i] =
663  send_adapt_child_num[i - 1] + send_adapt_child_num[i];
664  }
665 
666  send_adapt_child = (int *)calloc(send_adapt_child_num[pesize], sizeof(int));
667  send_index_child = (int *)calloc(send_elem_num[pesize], sizeof(int));
668  if ((send_adapt_child == NULL) || (send_index_child == NULL))
669  HECMW_dlb_memory_exit("send_adapt_child, send_index_child");
670 
671  for (i = 0; i < pesize + 1; i++) send_tmp_num[i] = 0;
672  tmp_int = -1;
673  for (i = 0; i < send_elem_num[pesize]; i++) {
674  for (j = mesh->adapt_children_index[send_elem[i]];
675  j < mesh->adapt_children_index[send_elem[i] + 1]; j++) {
676  tmp_int++;
677  if (mesh->adapt_children_item[j * 2 + 1] < 0)
678  send_adapt_child[tmp_int] = -1;
679  else if (mesh->adapt_children_item[j * 2 + 1] != mynode)
680  send_tmp_num[mesh->adapt_children_item[j * 2 + 1] + 1]++;
681  else if (new_elem[inter_elem[mesh->adapt_children_item[j * 2] - 1]] != -1)
682  send_adapt_child[tmp_int] =
683  new_elem[inter_elem[mesh->adapt_children_item[j * 2] - 1]];
684  else if (new_elem[inter_elem[mesh->adapt_children_item[j * 2] - 1]] ==
685  -1) {
686  fprintf(stderr, "There is something wrong with child information\n");
687  fprintf(stderr, "i=%d, send_elem[i]=%d child is %d PE=%d\n", i,
688  send_elem[i], mesh->adapt_children_item[j * 2] - 1,
689  mesh->adapt_children_item[j * 2 + 1]);
690  HECMW_dlb_print_exit("Error: in finding children inf");
691  }
692  }
693  }
694  for (i = 1; i < pesize + 1; i++) {
695  send_tmp_num[i] = send_tmp_num[i - 1] + send_tmp_num[i];
696  }
697 
698  stack_whole_send_recv(pesize, send_tmp_num, recv_tmp_num, mesh->HECMW_COMM,
699  mynode);
700  for (i = 0; i < pesize; i++) count_num[i] = 0;
701  if (send_tmp_num[pesize] > 0) {
702  send_tmp = (int *)calloc(send_tmp_num[pesize], sizeof(int));
703  if (send_tmp == NULL) HECMW_dlb_memory_exit("send_tmp");
704  }
705  if (recv_tmp_num[pesize] > 0) {
706  recv_tmp = (int *)calloc(recv_tmp_num[pesize], sizeof(int));
707  if (recv_tmp == NULL) HECMW_dlb_memory_exit("recv_tmp");
708  }
709  tmp_int = -1;
710  for (i = 0; i < send_elem_num[pesize]; i++) {
711  for (j = mesh->adapt_children_index[send_elem[i]];
712  j < mesh->adapt_children_index[send_elem[i] + 1]; j++) {
713  tmp_int++;
714  if (mesh->adapt_children_item[j * 2 + 1] < 0)
715  send_adapt_child[tmp_int] = -1;
716  else if (mesh->adapt_children_item[j * 2 + 1] != mynode) {
717  send_tmp[send_tmp_num[mesh->adapt_children_item[j * 2 + 1]] +
718  count_num[mesh->adapt_children_item[j * 2 + 1]]] =
719  inter_elem[mesh->adapt_children_item[j * 2] - 1];
720  count_num[mesh->adapt_children_item[j * 2 + 1]]++;
721  }
722  }
723  }
724  if (send_tmp_num[pesize] > 0)
725  int2_whole_send_recv(send_tmp_num[pesize], recv_tmp_num[pesize], pesize,
726  recv_tmp_num, send_tmp_num, send_tmp, recv_tmp,
727  mesh->HECMW_COMM, mynode);
728  for (i = 0; i < recv_tmp_num[pesize]; i++) {
729  if (new_elem[recv_tmp[i]] == -1)
730  HECMW_dlb_print_exit("There is something wrong in children inf");
731  else
732  recv_tmp[i] = new_elem[recv_tmp[i]];
733  }
734  if (send_tmp_num[pesize] > 0)
735  int2_whole_send_recv(recv_tmp_num[pesize], send_tmp_num[pesize], pesize,
736  send_tmp_num, recv_tmp_num, recv_tmp, send_tmp,
737  mesh->HECMW_COMM, mynode);
738 
739  for (i = 0; i < pesize; i++) count_num[i] = 0;
740  tmp_int = -1;
741  for (i = 0; i < send_elem_num[pesize]; i++) {
742  send_index_child[i] = mesh->adapt_children_index[send_elem[i] + 1] -
743  mesh->adapt_children_index[send_elem[i]];
744  for (j = mesh->adapt_children_index[send_elem[i]];
745  j < mesh->adapt_children_index[send_elem[i] + 1]; j++) {
746  tmp_int++;
747  if (mesh->adapt_children_item[j * 2 + 1] < 0) {
748  send_adapt_child[tmp_int] = -1;
749  } else if (mesh->adapt_children_item[j * 2 + 1] != mynode) {
750  send_adapt_child[tmp_int] =
751  send_tmp[send_tmp_num[mesh->adapt_children_item[j * 2 + 1]] +
752  count_num[mesh->adapt_children_item[j * 2 + 1]]];
753  count_num[mesh->adapt_children_item[j * 2 + 1]]++;
754  } else if (new_elem[inter_elem[mesh->adapt_children_item[j * 2] - 1]] !=
755  -1)
756  send_adapt_child[tmp_int] =
757  new_elem[inter_elem[mesh->adapt_children_item[j * 2] - 1]];
758  else if (new_elem[inter_elem[mesh->adapt_children_item[j * 2] - 1]] ==
759  -1) {
760  fprintf(stderr, "There is something wrong with child information\n");
761  fprintf(stderr, "i=%d, send_elem[i]=%d child is %d PE=%d\n", i,
762  send_elem[i], mesh->adapt_children_item[j * 2] - 1,
763  mesh->adapt_children_item[j * 2 + 1]);
764  HECMW_dlb_print_exit("There is something wrong in children inf");
765  }
766  }
767  }
768 
769  recv_adapt_child_num = (int *)calloc(pesize + 1, sizeof(int));
770  if (recv_adapt_child_num == NULL)
771  HECMW_dlb_memory_exit("recv_adapt_child_num");
772 
773  stack_whole_send_recv(pesize, send_adapt_child_num, recv_adapt_child_num,
774  mesh->HECMW_COMM, mynode);
775  recv_adapt_child = (int *)calloc(recv_adapt_child_num[pesize], sizeof(int));
776  if (recv_adapt_child == NULL) HECMW_dlb_memory_exit("recv_adapt_child");
778  send_adapt_child_num[pesize], recv_adapt_child_num[pesize], pesize,
779  recv_adapt_child_num, send_adapt_child_num, send_adapt_child,
780  recv_adapt_child, mesh->HECMW_COMM, mynode);
781 
782  int2_whole_send_recv(send_elem_num[pesize], recv_elem_num[pesize], pesize,
783  recv_elem_num, send_elem_num, send_index_child,
784  tmp_int_recv, mesh->HECMW_COMM, mynode);
785 
786  /*--------first find child number in send_parent --------------- */
787 
788  send2_adapt_child_num = (int *)calloc(pesize + 1, sizeof(int));
789  for (j = 0; j < pesize + 1; j++) send2_adapt_child_num[j] = 0;
790  for (j = 0; j < pesize; j++) {
791  for (i = send_parent_num[j]; i < send_parent_num[j + 1]; i++)
792  send2_adapt_child_num[j + 1] +=
793  mesh->adapt_children_index[send_parent[i] + 1] -
794  mesh->adapt_children_index[send_parent[i]];
795  }
796  for (i = 1; i < pesize + 1; i++) {
797  send2_adapt_child_num[i] =
798  send2_adapt_child_num[i - 1] + send2_adapt_child_num[i];
799  }
800  recv2_adapt_child_num = (int *)calloc(pesize + 1, sizeof(int));
801  if (recv2_adapt_child_num == NULL)
802  HECMW_dlb_memory_exit("recv2_adapt_child_num");
803 
804  stack_whole_send_recv(pesize, send2_adapt_child_num, recv2_adapt_child_num,
805  mesh->HECMW_COMM, mynode);
806  new_mesh->adapt_children_item = (int *)calloc(
807  2 * (recv_adapt_child_num[pesize] + recv2_adapt_child_num[pesize]),
808  sizeof(int));
810  HECMW_dlb_memory_exit("new_mesh: adapt_children_item");
812  (int *)calloc(new_mesh->n_elem + 1, sizeof(int));
814  HECMW_dlb_memory_exit("new_mesh: adapt_children_index");
815  /*
816  new_l_child=recv_adapt_child_num[pesize]+recv2_adapt_child_num[pesize];
817  */
818 
819  for (i = 0; i < recv_adapt_child_num[pesize]; i++) {
820  if (recv_adapt_child[i] == -1) {
821  new_mesh->adapt_children_item[i * 2] = 0;
822  new_mesh->adapt_children_item[i * 2 + 1] = -1;
823  } else {
824  new_mesh->adapt_children_item[i * 2] = (recv_adapt_child[i] % t_elem) + 1;
825  new_mesh->adapt_children_item[i * 2 + 1] = recv_adapt_child[i] / t_elem;
826  }
827  }
829  for (i = 0; i < recv_elem_num[pesize]; i++)
831  new_mesh->adapt_children_index[i] + tmp_int_recv[i];
832  if (send_tmp_num[pesize] > 0) free(send_tmp);
833  if (recv_tmp_num[pesize] > 0) free(recv_tmp);
834  if (send_elem_num[pesize] > 0) {
835  free(send_adapt_child);
836  free(recv_adapt_child);
837  free(send_index_child);
838  }
839 
840  send_adapt_child = (int *)calloc(send2_adapt_child_num[pesize], sizeof(int));
841  send_index_child = (int *)calloc(send_parent_num[pesize], sizeof(int));
842  if ((send_adapt_child == NULL) || (send_index_child == NULL))
843  HECMW_dlb_memory_exit("send_adapt_child, send_index_child");
844 
845  for (i = 0; i < pesize + 1; i++) send_tmp_num[i] = 0;
846  tmp_int = -1;
847  for (i = 0; i < send_parent_num[pesize]; i++) {
848  for (j = mesh->adapt_children_index[send_parent[i]];
849  j < mesh->adapt_children_index[send_parent[i] + 1]; j++) {
850  tmp_int++;
851  if (mesh->adapt_children_item[j * 2 + 1] < 0)
852  send_adapt_child[tmp_int] = -1;
853  else if (mesh->adapt_children_item[j * 2 + 1] != mynode)
854  send_tmp_num[mesh->adapt_children_item[j * 2 + 1] + 1]++;
855  else if (new_elem[inter_elem[mesh->adapt_children_item[j * 2] - 1]] != -1)
856  send_adapt_child[tmp_int] =
857  new_elem[inter_elem[mesh->adapt_children_item[j * 2] - 1]];
858  else if (new_elem[inter_elem[mesh->adapt_children_item[j * 2] - 1]] ==
859  -1) {
860  fprintf(stderr, "There is something wrong with child information\n");
861  fprintf(stderr, "in PE %d i=%d, send_parent[i]=%d child is %d PE=%d\n",
862  mynode, i, send_parent[i], mesh->adapt_children_item[j * 2] - 1,
863  mesh->adapt_children_item[j * 2 + 1]);
864  HECMW_dlb_print_exit("Error in finding children inf");
865  }
866  }
867  }
868  for (i = 1; i < pesize + 1; i++) {
869  send_tmp_num[i] = send_tmp_num[i - 1] + send_tmp_num[i];
870  }
871 
872  stack_whole_send_recv(pesize, send_tmp_num, recv_tmp_num, mesh->HECMW_COMM,
873  mynode);
874  for (i = 0; i < pesize; i++) count_num[i] = 0;
875  if (send_tmp_num[pesize] > 0) {
876  send_tmp = (int *)calloc(send_tmp_num[pesize], sizeof(int));
877  if (send_tmp == NULL) HECMW_dlb_memory_exit("send_tmp");
878  }
879  if (recv_tmp_num[pesize] > 0) {
880  recv_tmp = (int *)calloc(recv_tmp_num[pesize], sizeof(int));
881  if (recv_tmp == NULL) HECMW_dlb_memory_exit("recv_tmp");
882  }
883  tmp_int = -1;
884  for (i = 0; i < send_parent_num[pesize]; i++) {
885  for (j = mesh->adapt_children_index[send_parent[i]];
886  j < mesh->adapt_children_index[send_parent[i] + 1]; j++) {
887  tmp_int++;
888  if (mesh->adapt_children_item[j * 2 + 1] < 0)
889  send_adapt_child[tmp_int] = -1;
890  else if (mesh->adapt_children_item[j * 2 + 1] != mynode) {
891  send_tmp[send_tmp_num[mesh->adapt_children_item[j * 2 + 1]] +
892  count_num[mesh->adapt_children_item[j * 2 + 1]]] =
893  inter_elem[mesh->adapt_children_item[j * 2] - 1];
894  count_num[mesh->adapt_children_item[j * 2 + 1]]++;
895  }
896  }
897  }
898  if (send_tmp_num[pesize] > 0)
899  int2_whole_send_recv(send_tmp_num[pesize], recv_tmp_num[pesize], pesize,
900  recv_tmp_num, send_tmp_num, send_tmp, recv_tmp,
901  mesh->HECMW_COMM, mynode);
902  for (i = 0; i < recv_tmp_num[pesize]; i++) {
903  if (new_elem[recv_tmp[i]] == -1)
904  HECMW_dlb_print_exit("There is something wrong in parents' children inf");
905  else
906  recv_tmp[i] = new_elem[recv_tmp[i]];
907  }
908  if (send_tmp_num[pesize] > 0)
909  int2_whole_send_recv(recv_tmp_num[pesize], send_tmp_num[pesize], pesize,
910  send_tmp_num, recv_tmp_num, recv_tmp, send_tmp,
911  mesh->HECMW_COMM, mynode);
912 
913  for (i = 0; i < pesize; i++) count_num[i] = 0;
914  tmp_int = -1;
915  for (i = 0; i < send_parent_num[pesize]; i++) {
916  send_index_child[i] = mesh->adapt_children_index[send_parent[i] + 1] -
917  mesh->adapt_children_index[send_parent[i]];
918  for (j = mesh->adapt_children_index[send_parent[i]];
919  j < mesh->adapt_children_index[send_parent[i] + 1]; j++) {
920  tmp_int++;
921  if (mesh->adapt_children_item[j * 2 + 1] < 0) {
922  send_adapt_child[tmp_int] = -1;
923  } else if (mesh->adapt_children_item[j * 2 + 1] != mynode) {
924  send_adapt_child[tmp_int] =
925  send_tmp[send_tmp_num[mesh->adapt_children_item[j * 2 + 1]] +
926  count_num[mesh->adapt_children_item[j * 2 + 1]]];
927  count_num[mesh->adapt_children_item[j * 2 + 1]]++;
928  } else if (new_elem[inter_elem[mesh->adapt_children_item[j * 2] - 1]] !=
929  -1)
930  send_adapt_child[tmp_int] =
931  new_elem[inter_elem[mesh->adapt_children_item[j * 2] - 1]];
932  else if (new_elem[inter_elem[mesh->adapt_children_item[j * 2] - 1]] ==
933  -1) {
934  fprintf(stderr, "There is something wrong with child information\n");
935  fprintf(stderr, "i=%d, send_elem[i]=%d child is %d PE=%d\n", i,
936  send_parent[i], mesh->adapt_children_item[j * 2] - 1,
937  mesh->adapt_children_item[j * 2 + 1]);
938  HECMW_dlb_print_exit("ERROR in finding parent elements' children inf");
939  }
940  }
941  }
942 
943  recv_adapt_child = (int *)calloc(recv2_adapt_child_num[pesize], sizeof(int));
944  if (recv_adapt_child == NULL) HECMW_dlb_memory_exit("recv_adapt_child");
946  send2_adapt_child_num[pesize], recv2_adapt_child_num[pesize], pesize,
947  recv2_adapt_child_num, send2_adapt_child_num, send_adapt_child,
948  recv_adapt_child, mesh->HECMW_COMM, mynode);
949 
950  int2_whole_send_recv(send_parent_num[pesize], recv_parent_num[pesize], pesize,
951  recv_parent_num, send_parent_num, send_index_child,
952  tmp2_int_recv, mesh->HECMW_COMM, mynode);
953 
954  for (i = 0; i < recv2_adapt_child_num[pesize]; i++) {
955  if (recv_adapt_child[i] == -1) {
956  new_mesh->adapt_children_item[(i + recv_adapt_child_num[pesize]) * 2] = 0;
957  new_mesh
958  ->adapt_children_item[(i + recv_adapt_child_num[pesize]) * 2 + 1] =
959  -1;
960  } else {
961  new_mesh->adapt_children_item[(i + recv_adapt_child_num[pesize]) * 2] =
962  (recv_adapt_child[i] % t_elem) + 1;
963  new_mesh
964  ->adapt_children_item[(i + recv_adapt_child_num[pesize]) * 2 + 1] =
965  recv_adapt_child[i] / t_elem;
966  }
967  }
968  new_mesh->adapt_children_index[recv_elem_num[pesize] + 1] =
969  new_mesh->adapt_children_index[recv_elem_num[pesize]] + tmp2_int_recv[0];
970  for (i = 1; i < recv_parent_num[pesize]; i++)
971  new_mesh->adapt_children_index[i + 1 + recv_elem_num[pesize]] =
972  new_mesh->adapt_children_index[i + recv_elem_num[pesize]] +
973  tmp2_int_recv[i];
974  if (send_tmp_num[pesize] > 0) free(send_tmp);
975  if (recv_tmp_num[pesize] > 0) free(recv_tmp);
976  if (send_parent_num[pesize] > 0) {
977  free(send_adapt_child);
978  free(recv_adapt_child);
979  free(send_index_child);
980  }
981 
982  /*---------------send elem_id --------------------------------*/
983 
984  if (new_mesh->n_elem > 0) {
985  new_mesh->elem_ID = (int *)calloc(2 * new_mesh->n_elem, sizeof(int));
986  if (new_mesh->elem_ID == NULL) HECMW_dlb_memory_exit("new_mesh: elem_id");
987  }
988 
989  tmp_int_send = (int *)calloc(send_elem_num[pesize] + 1, sizeof(int));
990 
991  if (tmp_int_send == NULL) HECMW_dlb_memory_exit("tmp_int_send ");
992  if (send_elem_num[pesize] > 0) {
993  for (i = 0; i < send_elem_num[pesize]; i++) tmp_int_send[i] = send_inter[i];
994  }
995  tmp_int_recv = (int *)calloc(recv_elem_num[pesize] + 1, sizeof(int));
996  if (tmp_int_recv == NULL) HECMW_dlb_memory_exit("tmp_int_recv");
997 
998  int2_whole_send_recv(send_elem_num[pesize], recv_elem_num[pesize], pesize,
999  recv_elem_num, send_elem_num, tmp_int_send, tmp_int_recv,
1000  mesh->HECMW_COMM, mynode);
1001  if (recv_elem_num[pesize] > 0) {
1002  for (i = 0; i < recv_elem_num[pesize]; i++) {
1003  new_mesh->elem_ID[i * 2] = (tmp_int_recv[i] % t_elem) + 1;
1004  new_mesh->elem_ID[i * 2 + 1] = tmp_int_recv[i] / t_elem;
1005  }
1006  }
1007  tmp2_int_send = (int *)calloc(send_parent_num[pesize], sizeof(int));
1008  if (tmp2_int_send == NULL) HECMW_dlb_memory_exit("tmp2_int_send");
1009 
1010  for (i = 0; i < send_parent_num[pesize]; i++)
1011  tmp2_int_send[i] = send_parent_inter[i];
1012  int2_whole_send_recv(send_parent_num[pesize], recv_parent_num[pesize], pesize,
1013  recv_parent_num, send_parent_num, tmp2_int_send,
1014  tmp2_int_recv, mesh->HECMW_COMM, mynode);
1015  for (i = 0; i < recv_parent_num[pesize]; i++) {
1016  new_mesh->elem_ID[(i + recv_elem_num[pesize]) * 2] =
1017  (tmp2_int_recv[i] % t_elem) + 1;
1018  new_mesh->elem_ID[(i + recv_elem_num[pesize]) * 2 + 1] =
1019  tmp2_int_recv[i] / t_elem;
1020  }
1021 
1022  if (mynode == 0) fprintf(stderr, "Finish sending elem_id\n");
1023  if (mesh->n_elem > 0) {
1024  free(new_elem);
1025  free(send_inter);
1026  free(send_parent_inter);
1027  }
1028  /* ------ generate elem_internal_list and global_elem_id -------------- */
1029  new_mesh->ne_internal = recv_inter_num[pesize] + recv_parent_num[pesize];
1031  (int *)calloc(new_mesh->ne_internal, sizeof(int));
1033  HECMW_dlb_memory_exit("new_mesh: elem_internal_list");
1034  new_nelem_dist = (int *)calloc(pesize + 1, sizeof(int));
1035  if (new_nelem_dist == NULL) HECMW_dlb_memory_exit("new_nelem_dist");
1036  if (mynode == 0) {
1037  new_nelem_dist[0] = 0;
1038  new_nelem_dist[1] = new_mesh->ne_internal;
1039  tmp_sum = new_mesh->ne_internal;
1040  for (i = 1; i < pesize; i++) {
1041  HECMW_Recv(&tmp_int, 1, HECMW_INT, i, HECMW_ANY_TAG, mesh->HECMW_COMM,
1042  &stat);
1043  tmp_sum += tmp_int;
1044  new_nelem_dist[i + 1] = tmp_sum;
1045  }
1046  for (i = 1; i < pesize; i++)
1047  HECMW_Send(new_nelem_dist, pesize + 1, HECMW_INT, i, 0, mesh->HECMW_COMM);
1048  } else {
1050  HECMW_Recv(new_nelem_dist, pesize + 1, HECMW_INT, 0, HECMW_ANY_TAG,
1051  mesh->HECMW_COMM, &stat);
1052  }
1053  new_mesh->global_elem_ID = (int *)calloc(new_mesh->n_elem, sizeof(int));
1054  if (new_mesh->global_elem_ID == NULL)
1055  HECMW_dlb_memory_exit("new_mesh: global_elem_ID");
1056 
1057  tmp_int = 0;
1058  for (i = 0; i < new_mesh->n_elem; i++) {
1059  new_mesh->global_elem_ID[i] =
1060  new_nelem_dist[new_mesh->elem_ID[i * 2 + 1]] + new_mesh->elem_ID[i * 2];
1061  if (new_mesh->elem_ID[i * 2 + 1] == mynode) {
1062  new_mesh->elem_internal_list[tmp_int] = i + 1;
1063  tmp_int++;
1064  }
1065  }
1066 
1067  /*---------------send elem_type --------------------------------*/
1068  if (new_mesh->n_elem > 0) {
1069  new_mesh->elem_type = (int *)calloc(new_mesh->n_elem, sizeof(int));
1070  if (new_mesh->elem_type == NULL)
1071  HECMW_dlb_memory_exit("new_mesh: elem_type");
1072  }
1073  if (send_elem_num[pesize] > 0) {
1074  for (i = 0; i < send_elem_num[pesize]; i++)
1075  tmp_int_send[i] = mesh->elem_type[send_elem[i]];
1076  }
1077  int2_whole_send_recv(send_elem_num[pesize], recv_elem_num[pesize], pesize,
1078  recv_elem_num, send_elem_num, tmp_int_send, tmp_int_recv,
1079  mesh->HECMW_COMM, mynode);
1080  for (i = 0; i < recv_elem_num[pesize]; i++)
1081  new_mesh->elem_type[i] = tmp_int_recv[i];
1082  for (i = 0; i < send_parent_num[pesize]; i++)
1083  tmp2_int_send[i] = mesh->elem_type[send_parent[i]];
1084  int2_whole_send_recv(send_parent_num[pesize], recv_parent_num[pesize], pesize,
1085  recv_parent_num, send_parent_num, tmp2_int_send,
1086  tmp2_int_recv, mesh->HECMW_COMM, mynode);
1087  for (i = 0; i < recv_parent_num[pesize]; i++)
1088  new_mesh->elem_type[i + recv_elem_num[pesize]] = tmp2_int_recv[i];
1089  if (mynode == 0) fprintf(stderr, "Finish sending elem_type\n");
1090 
1091  if (new_mesh->n_elem > 0) {
1092  new_mesh->section_ID = (int *)calloc(new_mesh->n_elem, sizeof(int));
1093  if (new_mesh->section_ID == NULL)
1094  HECMW_dlb_memory_exit("new_mesh: section_ID");
1095  }
1096  if (send_elem_num[pesize] > 0) {
1097  for (i = 0; i < send_elem_num[pesize]; i++)
1098  tmp_int_send[i] = mesh->section_ID[send_elem[i]];
1099  }
1100  int2_whole_send_recv(send_elem_num[pesize], recv_elem_num[pesize], pesize,
1101  recv_elem_num, send_elem_num, tmp_int_send, tmp_int_recv,
1102  mesh->HECMW_COMM, mynode);
1103  for (i = 0; i < recv_elem_num[pesize]; i++)
1104  new_mesh->section_ID[i] = tmp_int_recv[i];
1105  for (i = 0; i < send_parent_num[pesize]; i++)
1106  tmp2_int_send[i] = mesh->section_ID[send_parent[i]];
1107  int2_whole_send_recv(send_parent_num[pesize], recv_parent_num[pesize], pesize,
1108  recv_parent_num, send_parent_num, tmp2_int_send,
1109  tmp2_int_recv, mesh->HECMW_COMM, mynode);
1110  for (i = 0; i < recv_parent_num[pesize]; i++)
1111  new_mesh->section_ID[i + recv_elem_num[pesize]] = tmp2_int_recv[i];
1112  if (mynode == 0) fprintf(stderr, "Finish sending section_ID\n");
1113 
1114  /* ---------- send material inf. ---------------- */
1115  if (new_mesh->n_elem > 0) {
1116  if (mesh->elem_mat_int_val != NULL) {
1118  (double *)calloc(new_mesh->n_elem, sizeof(double));
1119  if (new_mesh->elem_mat_int_val == NULL)
1120  HECMW_dlb_memory_exit("new_mesh: elem_mat_int_val");
1121  tmp_send_d = (double *)calloc(send_elem_num[pesize] + 1, sizeof(double));
1122  tmp_recv_d = (double *)calloc(recv_elem_num[pesize] + 1, sizeof(double));
1123  if (tmp_send_d == NULL) HECMW_dlb_memory_exit("tmp_send_d");
1124  if (send_elem_num[pesize] > 0) {
1125  for (i = 0; i < send_elem_num[pesize]; i++)
1126  tmp_send_d[i] = mesh->elem_mat_int_val[send_elem[i]];
1127  }
1128  double2_whole_send_recv(send_elem_num[pesize], recv_elem_num[pesize],
1129  pesize, recv_elem_num, send_elem_num, tmp_send_d,
1130  tmp_recv_d, mesh->HECMW_COMM, mynode);
1131  for (i = 0; i < recv_elem_num[pesize]; i++)
1132  new_mesh->elem_mat_int_val[i] = tmp_recv_d[i];
1133  free(tmp_send_d);
1134  free(tmp_recv_d);
1135  tmp2_send_d =
1136  (double *)calloc(send_parent_num[pesize] + 1, sizeof(double));
1137  tmp2_recv_d =
1138  (double *)calloc(recv_parent_num[pesize] + 1, sizeof(double));
1139  if (tmp2_send_d == NULL) HECMW_dlb_memory_exit("tmp2_send_d");
1140 
1141  for (i = 0; i < send_parent_num[pesize]; i++)
1142  tmp2_send_d[i] = mesh->elem_mat_int_val[send_parent[i]];
1143  double2_whole_send_recv(send_parent_num[pesize], recv_parent_num[pesize],
1144  pesize, recv_parent_num, send_parent_num,
1145  tmp2_send_d, tmp2_recv_d, mesh->HECMW_COMM,
1146  mynode);
1147  for (i = 0; i < recv_parent_num[pesize]; i++)
1148  new_mesh->elem_mat_int_val[i + recv_elem_num[pesize]] = tmp2_recv_d[i];
1149  free(tmp2_send_d);
1150  free(tmp2_recv_d);
1151  }
1152  }
1153  if (new_mesh->n_elem > 0) {
1155  (int *)calloc(new_mesh->n_elem + 1, sizeof(int));
1156  new_mesh->elem_mat_ID_item = (int *)calloc(new_mesh->n_elem, sizeof(int));
1157  if ((new_mesh->elem_mat_ID_index == NULL) ||
1159  HECMW_dlb_memory_exit("new_mesh: elem_mat_ID_index, elem_mat_ID_item");
1160  new_mesh->elem_mat_ID_index[0] = 0;
1161  for (i = 0; i < new_mesh->n_elem; i++) {
1162  new_mesh->elem_mat_ID_item[i] = 1;
1163  new_mesh->elem_mat_ID_index[i + 1] = i + 1;
1164  }
1165  }
1166 
1167  if (mynode == 0) fprintf(stderr, "Finish sending elem_material inf\n");
1168 
1169  /*---------------send adaptation_level --------------------------------*/
1172  new_mesh->adapt_level = (int *)calloc(new_mesh->n_elem, sizeof(int));
1173  if (new_mesh->adapt_level == NULL)
1174  HECMW_dlb_memory_exit("new_mesh: adapt_level");
1175 
1176  for (i = 0; i < send_elem_num[pesize]; i++)
1177  tmp_int_send[i] = mesh->adapt_level[send_elem[i]];
1178  int2_whole_send_recv(send_elem_num[pesize], recv_elem_num[pesize], pesize,
1179  recv_elem_num, send_elem_num, tmp_int_send, tmp_int_recv,
1180  mesh->HECMW_COMM, mynode);
1181  for (i = 0; i < recv_elem_num[pesize]; i++)
1182  new_mesh->adapt_level[i] = tmp_int_recv[i];
1183 
1184  for (i = 0; i < send_parent_num[pesize]; i++)
1185  tmp2_int_send[i] = mesh->adapt_level[send_parent[i]];
1186  int2_whole_send_recv(send_parent_num[pesize], recv_parent_num[pesize], pesize,
1187  recv_parent_num, send_parent_num, tmp2_int_send,
1188  tmp2_int_recv, mesh->HECMW_COMM, mynode);
1189  for (i = 0; i < recv_parent_num[pesize]; i++)
1190  new_mesh->adapt_level[i + recv_elem_num[pesize]] = tmp2_int_recv[i];
1191  if (mynode == 0) fprintf(stderr, "Finish sending adaptation_level\n");
1192  /*---------------send adaptation_type --------------------------------*/
1193  new_mesh->adapt_type = (int *)calloc(new_mesh->n_elem, sizeof(int));
1194  if (new_mesh->adapt_type == NULL)
1195  HECMW_dlb_memory_exit("new_mesh: adapt_type");
1196 
1197  for (i = 0; i < send_elem_num[pesize]; i++)
1198  tmp_int_send[i] = mesh->adapt_type[send_elem[i]];
1199  int2_whole_send_recv(send_elem_num[pesize], recv_elem_num[pesize], pesize,
1200  recv_elem_num, send_elem_num, tmp_int_send, tmp_int_recv,
1201  mesh->HECMW_COMM, mynode);
1202  for (i = 0; i < recv_elem_num[pesize]; i++)
1203  new_mesh->adapt_type[i] = tmp_int_recv[i];
1204 
1205  for (i = 0; i < send_parent_num[pesize]; i++)
1206  tmp2_int_send[i] = mesh->adapt_type[send_parent[i]];
1207  int2_whole_send_recv(send_parent_num[pesize], recv_parent_num[pesize], pesize,
1208  recv_parent_num, send_parent_num, tmp2_int_send,
1209  tmp2_int_recv, mesh->HECMW_COMM, mynode);
1210  for (i = 0; i < recv_parent_num[pesize]; i++)
1211  new_mesh->adapt_type[i + recv_elem_num[pesize]] = tmp2_int_recv[i];
1212  if (mynode == 0) fprintf(stderr, "Finish sending adaptation_type\n");
1213 
1214  /* send wheniwasrefined elem ------------- */
1216  (int *)calloc(new_mesh->n_elem, sizeof(int));
1218  HECMW_dlb_memory_exit("new_mesh: when_i_was_refined_elem");
1219 
1220  for (i = 0; i < send_elem_num[pesize]; i++)
1221  tmp_int_send[i] = mesh->when_i_was_refined_elem[send_elem[i]];
1222  int2_whole_send_recv(send_elem_num[pesize], recv_elem_num[pesize], pesize,
1223  recv_elem_num, send_elem_num, tmp_int_send, tmp_int_recv,
1224  mesh->HECMW_COMM, mynode);
1225  for (i = 0; i < recv_elem_num[pesize]; i++)
1226  new_mesh->when_i_was_refined_elem[i] = tmp_int_recv[i];
1227 
1228  for (i = 0; i < send_parent_num[pesize]; i++)
1229  tmp2_int_send[i] = mesh->when_i_was_refined_elem[send_parent[i]];
1230  int2_whole_send_recv(send_parent_num[pesize], recv_parent_num[pesize], pesize,
1231  recv_parent_num, send_parent_num, tmp2_int_send,
1232  tmp2_int_recv, mesh->HECMW_COMM, mynode);
1233  for (i = 0; i < recv_parent_num[pesize]; i++)
1234  new_mesh->when_i_was_refined_elem[i + recv_elem_num[pesize]] =
1235  tmp2_int_recv[i];
1236  if (mynode == 0) fprintf(stderr, "Finish sending when_i_was_refined_elem\n");
1237 
1238  /* ------- send elem_index -------- */
1239  if (new_mesh->n_elem > 0) {
1241  (int *)calloc(new_mesh->n_elem + 1, sizeof(int));
1242  if (new_mesh->elem_node_index == NULL)
1243  HECMW_dlb_memory_exit("new_mesh: elem_node_index");
1244  }
1245  for (i = 0; i < send_elem_num[pesize]; i++)
1246  tmp_int_send[i] = mesh->elem_node_index[send_elem[i] + 1] -
1247  mesh->elem_node_index[send_elem[i]];
1248  int2_whole_send_recv(send_elem_num[pesize], recv_elem_num[pesize], pesize,
1249  recv_elem_num, send_elem_num, tmp_int_send, tmp_int_recv,
1250  mesh->HECMW_COMM, mynode);
1251  for (i = 0; i < send_parent_num[pesize]; i++)
1252  tmp2_int_send[i] = mesh->elem_node_index[send_parent[i] + 1] -
1253  mesh->elem_node_index[send_parent[i]];
1254  int2_whole_send_recv(send_parent_num[pesize], recv_parent_num[pesize], pesize,
1255  recv_parent_num, send_parent_num, tmp2_int_send,
1256  tmp2_int_recv, mesh->HECMW_COMM, mynode);
1257 
1258  if (new_mesh->n_elem > 0) {
1259  new_mesh->elem_node_index[0] = 0;
1260  for (i = 0; i < recv_elem_num[pesize]; i++)
1261  new_mesh->elem_node_index[i + 1] =
1262  new_mesh->elem_node_index[i] + tmp_int_recv[i];
1263  for (i = 0; i < recv_parent_num[pesize]; i++)
1264  new_mesh->elem_node_index[recv_elem_num[pesize] + 1 + i] =
1265  tmp2_int_recv[i] +
1266  new_mesh->elem_node_index[recv_elem_num[pesize] + i];
1267  }
1268 
1269  if (mynode == 0) fprintf(stderr, "Finish sending elem_node_index\n");
1270  /* -------------- send elem_node_item -------------------------------------
1271  */
1272  send_ptr_num = (int *)calloc(pesize + 1, sizeof(int));
1273  send_ptr_parent_num = (int *)calloc(pesize + 1, sizeof(int));
1274  recv_ptr_num = (int *)calloc(pesize + 1, sizeof(int));
1275  recv_ptr_parent_num = (int *)calloc(pesize + 1, sizeof(int));
1276  if ((send_ptr_num == NULL) || (send_ptr_parent_num == NULL) ||
1277  (recv_ptr_num == NULL) || (recv_ptr_parent_num == NULL))
1278  HECMW_dlb_memory_exit("send_recv_ptr_num, send_recv_parent_ptr_num");
1279  for (i = 0; i < pesize + 1; i++) {
1280  send_ptr_num[i] = 0;
1281  send_ptr_parent_num[i] = 0;
1282  }
1283  for (i = 1; i < pesize + 1; i++) {
1284  for (j = send_elem_num[i - 1]; j < send_elem_num[i]; j++)
1285  send_ptr_num[i] += tmp_int_send[j];
1286  }
1287  for (i = 1; i < pesize + 1; i++)
1288  send_ptr_num[i] = send_ptr_num[i - 1] + send_ptr_num[i];
1289  for (i = 1; i < pesize + 1; i++) {
1290  for (j = send_parent_num[i - 1]; j < send_parent_num[i]; j++)
1291  send_ptr_parent_num[i] += tmp2_int_send[j];
1292  }
1293  for (i = 1; i < pesize + 1; i++)
1294  send_ptr_parent_num[i] =
1295  send_ptr_parent_num[i - 1] + send_ptr_parent_num[i];
1296 
1297  stack_whole_send_recv(pesize, send_ptr_num, recv_ptr_num, mesh->HECMW_COMM,
1298  mynode);
1299  stack_whole_send_recv(pesize, send_ptr_parent_num, recv_ptr_parent_num,
1300  mesh->HECMW_COMM, mynode);
1301 
1302  new_mesh->node_group =
1303  (struct hecmwST_node_grp *)calloc(1, sizeof(struct hecmwST_node_grp));
1304  new_mesh->elem_group =
1305  (struct hecmwST_elem_grp *)calloc(1, sizeof(struct hecmwST_elem_grp));
1306  new_mesh->surf_group =
1307  (struct hecmwST_surf_grp *)calloc(1, sizeof(struct hecmwST_surf_grp));
1308 
1312  if (mesh->node_group->n_grp > 0) {
1314  (char **)calloc(mesh->node_group->n_grp, sizeof(char *));
1315  for (m = 0; m < mesh->node_group->n_grp; m++) {
1316  new_mesh->node_group->grp_name[m] = (char *)calloc(128, sizeof(char));
1317  if (new_mesh->node_group->grp_name[m] == NULL)
1318  HECMW_dlb_memory_exit("new_mesh: grp_name");
1319  strcpy(new_mesh->node_group->grp_name[m], mesh->node_group->grp_name[m]);
1320  }
1321  }
1322  if (mesh->elem_group->n_grp > 0) {
1324  (char **)calloc(mesh->elem_group->n_grp, sizeof(char *));
1325  for (m = 0; m < mesh->elem_group->n_grp; m++) {
1326  new_mesh->elem_group->grp_name[m] = (char *)calloc(128, sizeof(char));
1327  if (new_mesh->elem_group->grp_name[m] == NULL)
1328  HECMW_dlb_memory_exit("new_mesh: grp_name");
1329  strcpy(new_mesh->elem_group->grp_name[m], mesh->elem_group->grp_name[m]);
1330  }
1331  }
1332  if (mesh->surf_group->n_grp > 0) {
1334  (char **)calloc(mesh->surf_group->n_grp, sizeof(char *));
1335  for (m = 0; m < mesh->surf_group->n_grp; m++) {
1336  new_mesh->surf_group->grp_name[m] = (char *)calloc(128, sizeof(char));
1337  if (new_mesh->surf_group->grp_name[m] == NULL)
1338  HECMW_dlb_memory_exit("new_mesh: grp_name");
1339  strcpy(new_mesh->surf_group->grp_name[m], mesh->surf_group->grp_name[m]);
1340  }
1341  }
1342 
1343  if (new_mesh->elem_group->n_grp > 0) {
1344  tmp_grp =
1345  (Tmp_grp_inf *)calloc(new_mesh->elem_group->n_grp, sizeof(Tmp_grp_inf));
1346  if (tmp_grp == NULL) HECMW_dlb_memory_exit("tmp_grp");
1347  /*
1348  tmp_int_send=(int *)calloc(send_elem_num[pesize]+1,
1349  sizeof(int));
1350  if(tmp_int_send==NULL)
1351  HECMW_dlb_memory_exit("tmp_int_send ");
1352  tmp_int_recv=(int *)calloc(recv_elem_num[pesize]+1,
1353  sizeof(int));
1354  if(tmp_int_recv==NULL)
1355  HECMW_dlb_memory_exit("tmp_int_recv");
1356  tmp2_int_send=(int *)calloc(send_parent_num[pesize]+1,
1357  sizeof(int));
1358  if(tmp2_int_send==NULL)
1359  HECMW_dlb_memory_exit("tmp2_int_send ");
1360  tmp2_int_recv=(int *)calloc(recv_parent_num[pesize]+1,
1361  sizeof(int));
1362  if(tmp2_int_recv==NULL)
1363  HECMW_dlb_memory_exit("tmp2_int_recv");
1364  */
1365 
1366  for (m = 0; m < new_mesh->elem_group->n_grp; m++)
1367  tmp_grp[m].num_of_item = 0;
1368  for (m = 0; m < new_mesh->elem_group->n_grp; m++) {
1369  if (m == 0) {
1370  tmp_elem_grp = (int *)calloc(mesh->n_elem, sizeof(int));
1371  if (tmp_elem_grp == NULL) HECMW_dlb_memory_exit("tmp_elem_grp");
1372  }
1373  for (i = 0; i < mesh->n_elem; i++) tmp_elem_grp[i] = 0;
1374  if ((mesh->elem_group->grp_index[m + 1] -
1375  mesh->elem_group->grp_index[m]) > 0) {
1376  for (i = mesh->elem_group->grp_index[m];
1377  i < mesh->elem_group->grp_index[m + 1]; i++)
1378  tmp_elem_grp[mesh->elem_group->grp_item[i] - 1] = 1;
1379  }
1380  if (send_elem_num[pesize] > 0) {
1381  for (i = 0; i < send_elem_num[pesize]; i++)
1382  tmp_int_send[i] = tmp_elem_grp[send_elem[i]];
1383  int2_whole_send_recv(send_elem_num[pesize], recv_elem_num[pesize],
1384  pesize, recv_elem_num, send_elem_num, tmp_int_send,
1385  tmp_int_recv, mesh->HECMW_COMM, mynode);
1386  }
1387  if (send_parent_num[pesize] > 0) {
1388  for (i = 0; i < send_parent_num[pesize]; i++)
1389  tmp2_int_send[i] = tmp_elem_grp[send_parent[i]];
1390  int2_whole_send_recv(send_parent_num[pesize], recv_parent_num[pesize],
1391  pesize, recv_parent_num, send_parent_num,
1392  tmp2_int_send, tmp2_int_recv, mesh->HECMW_COMM,
1393  mynode);
1394  }
1395  num_grp_item = 0;
1396  for (i = 0; i < recv_elem_num[pesize]; i++) {
1397  if (tmp_int_recv[i] == 1) num_grp_item++;
1398  }
1399  for (i = 0; i < recv_parent_num[pesize]; i++) {
1400  if (tmp2_int_recv[i] == 1) num_grp_item++;
1401  }
1402  tmp_grp[m].num_of_item = num_grp_item;
1403  if (num_grp_item > 0) {
1404  tmp_grp[m].item = (int *)calloc(num_grp_item, sizeof(int));
1405  if (tmp_grp[m].item == NULL) HECMW_dlb_memory_exit("tmp_grp:item");
1406  tmp_int = 0;
1407  for (i = 0; i < recv_elem_num[pesize]; i++) {
1408  if (tmp_int_recv[i] == 1) {
1409  tmp_grp[m].item[tmp_int] = i + 1;
1410  tmp_int++;
1411  }
1412  }
1413  for (i = 0; i < recv_parent_num[pesize]; i++) {
1414  if (tmp2_int_recv[i] == 1) {
1415  tmp_grp[m].item[tmp_int] = recv_elem_num[pesize] + i + 1;
1416  tmp_int++;
1417  }
1418  }
1419  }
1420  }
1421 
1422  free(tmp_elem_grp);
1423  num_grp_item = 0;
1424  for (m = 0; m < new_mesh->elem_group->n_grp; m++)
1425  num_grp_item += tmp_grp[m].num_of_item;
1427  (int *)calloc(new_mesh->elem_group->n_grp + 1, sizeof(int));
1428  if (new_mesh->elem_group->grp_index == NULL)
1429  HECMW_dlb_memory_exit("new_mesh: elem_grp: grp_index");
1430  new_mesh->elem_group->grp_index[0] = 0;
1431  for (m = 0; m < new_mesh->elem_group->n_grp; m++)
1432  new_mesh->elem_group->grp_index[m + 1] =
1433  new_mesh->elem_group->grp_index[m] + tmp_grp[m].num_of_item;
1434  if (num_grp_item > 0) {
1435  new_mesh->elem_group->grp_item = (int *)calloc(num_grp_item, sizeof(int));
1436  if (new_mesh->elem_group->grp_item == NULL)
1437  HECMW_dlb_memory_exit("new_mesh: elem_grp: grp_item");
1438  tmp_int = 0;
1439  for (m = 0; m < new_mesh->elem_group->n_grp; m++) {
1440  for (i = 0; i < tmp_grp[m].num_of_item; i++) {
1441  new_mesh->elem_group->grp_item[tmp_int] = tmp_grp[m].item[i];
1442  tmp_int++;
1443  }
1444  }
1445  }
1446  for (m = 0; m < new_mesh->elem_group->n_grp; m++) {
1447  if (tmp_grp[m].num_of_item > 0) free(tmp_grp[m].item);
1448  }
1449  free(tmp_grp);
1450  if (mynode == 0) fprintf(stderr, "Finish generating new elem_grp inf.\n");
1451  }
1452 
1453  if (new_mesh->surf_group->n_grp > 0) {
1454  tmp_grp =
1455  (Tmp_grp_inf *)calloc(new_mesh->surf_group->n_grp, sizeof(Tmp_grp_inf));
1456  if (tmp_grp == NULL) HECMW_dlb_memory_exit("tmp_grp");
1457 
1458  tmp_surf_id = (int *)calloc(recv_elem_num[pesize] + 1, sizeof(int));
1459  if (tmp_surf_id == NULL) HECMW_dlb_memory_exit("tmp_surf_id");
1460  tmp2_surf_id = (int *)calloc(recv_parent_num[pesize] + 1, sizeof(int));
1461  if (tmp2_surf_id == NULL) HECMW_dlb_memory_exit("tmp2_surf_id");
1462 
1463  for (m = 0; m < new_mesh->surf_group->n_grp; m++)
1464  tmp_grp[m].num_of_item = 0;
1465  for (m = 0; m < new_mesh->surf_group->n_grp; m++) {
1466  if (m == 0) {
1467  tmp_elem_grp = (int *)calloc(mesh->n_elem, sizeof(int));
1468  if (tmp_elem_grp == NULL) HECMW_dlb_memory_exit("tmp_elem_grp");
1469  tmp_surf_grp = (int *)calloc(mesh->n_elem, sizeof(int));
1470  if (tmp_surf_grp == NULL) HECMW_dlb_memory_exit("tmp_surf_grp");
1471  }
1472  for (i = 0; i < mesh->n_elem; i++) tmp_elem_grp[i] = 0;
1473  for (i = 0; i < mesh->n_elem; i++) tmp_surf_grp[i] = -1;
1474  if ((mesh->surf_group->grp_index[m + 1] -
1475  mesh->surf_group->grp_index[m]) > 0) {
1476  for (i = mesh->surf_group->grp_index[m];
1477  i < mesh->surf_group->grp_index[m + 1]; i++) {
1478  tmp_elem_grp[mesh->surf_group->grp_item[i * 2] - 1] = 1;
1479  tmp_surf_grp[mesh->surf_group->grp_item[i * 2] - 1] =
1480  mesh->surf_group->grp_item[i * 2 + 1];
1481  }
1482  }
1483  if (send_elem_num[pesize] > 0) {
1484  for (i = 0; i < send_elem_num[pesize]; i++)
1485  tmp_int_send[i] = tmp_elem_grp[send_elem[i]];
1486  int2_whole_send_recv(send_elem_num[pesize], recv_elem_num[pesize],
1487  pesize, recv_elem_num, send_elem_num, tmp_int_send,
1488  tmp_int_recv, mesh->HECMW_COMM, mynode);
1489  for (i = 0; i < send_elem_num[pesize]; i++)
1490  tmp_int_send[i] = tmp_surf_grp[send_elem[i]];
1491  int2_whole_send_recv(send_elem_num[pesize], recv_elem_num[pesize],
1492  pesize, recv_elem_num, send_elem_num, tmp_int_send,
1493  tmp_surf_id, mesh->HECMW_COMM, mynode);
1494  }
1495  if (send_parent_num[pesize] > 0) {
1496  for (i = 0; i < send_parent_num[pesize]; i++)
1497  tmp2_int_send[i] = tmp_elem_grp[send_parent[i]];
1498  int2_whole_send_recv(send_parent_num[pesize], recv_parent_num[pesize],
1499  pesize, recv_parent_num, send_parent_num,
1500  tmp2_int_send, tmp2_int_recv, mesh->HECMW_COMM,
1501  mynode);
1502  for (i = 0; i < send_parent_num[pesize]; i++)
1503  tmp2_int_send[i] = tmp_surf_grp[send_parent[i]];
1504  int2_whole_send_recv(send_parent_num[pesize], recv_parent_num[pesize],
1505  pesize, recv_parent_num, send_parent_num,
1506  tmp2_int_send, tmp2_surf_id, mesh->HECMW_COMM,
1507  mynode);
1508  }
1509  num_grp_item = 0;
1510  for (i = 0; i < recv_elem_num[pesize]; i++) {
1511  if (tmp_int_recv[i] == 1) num_grp_item++;
1512  }
1513  for (i = 0; i < recv_parent_num[pesize]; i++) {
1514  if (tmp2_int_recv[i] == 1) num_grp_item++;
1515  }
1516  tmp_grp[m].num_of_item = num_grp_item;
1517  if (num_grp_item > 0) {
1518  tmp_grp[m].item = (int *)calloc(num_grp_item * 2, sizeof(int));
1519  if (tmp_grp[m].item == NULL) HECMW_dlb_memory_exit("tmp_grp:item");
1520  tmp_int = 0;
1521  for (i = 0; i < recv_elem_num[pesize]; i++) {
1522  if (tmp_int_recv[i] == 1) {
1523  tmp_grp[m].item[tmp_int * 2] = i + 1;
1524  tmp_grp[m].item[tmp_int * 2 + 1] = tmp_surf_id[i];
1525  tmp_int++;
1526  }
1527  }
1528  for (i = 0; i < recv_parent_num[pesize]; i++) {
1529  if (tmp2_int_recv[i] == 1) {
1530  tmp_grp[m].item[tmp_int * 2] = recv_elem_num[pesize] + i + 1;
1531  tmp_grp[m].item[tmp_int * 2 + 1] = tmp2_surf_id[i];
1532 
1533  tmp_int++;
1534  }
1535  }
1536  }
1537  }
1538 
1539  free(tmp_elem_grp);
1540  free(tmp_surf_grp);
1541  free(tmp2_surf_id);
1542  free(tmp_surf_id);
1543  num_grp_item = 0;
1544  for (m = 0; m < new_mesh->surf_group->n_grp; m++)
1545  num_grp_item += tmp_grp[m].num_of_item;
1547  (int *)calloc(new_mesh->surf_group->n_grp + 1, sizeof(int));
1548  if (new_mesh->surf_group->grp_index == NULL)
1549  HECMW_dlb_memory_exit("new_mesh: surf_grp: grp_index");
1550  new_mesh->surf_group->grp_index[0] = 0;
1551  for (m = 0; m < new_mesh->surf_group->n_grp; m++)
1552  new_mesh->surf_group->grp_index[m + 1] =
1553  new_mesh->surf_group->grp_index[m] + tmp_grp[m].num_of_item;
1554  if (num_grp_item > 0) {
1556  (int *)calloc(num_grp_item * 2, sizeof(int));
1557  if (new_mesh->surf_group->grp_item == NULL)
1558  HECMW_dlb_memory_exit("new_mesh: surf_grp: grp_item");
1559  tmp_int = 0;
1560  for (m = 0; m < new_mesh->surf_group->n_grp; m++) {
1561  for (i = 0; i < tmp_grp[m].num_of_item; i++) {
1562  new_mesh->surf_group->grp_item[tmp_int * 2] = tmp_grp[m].item[i * 2];
1563  new_mesh->surf_group->grp_item[tmp_int * 2 + 1] =
1564  tmp_grp[m].item[i * 2 + 1];
1565  tmp_int++;
1566  }
1567  }
1568  }
1569  for (m = 0; m < new_mesh->surf_group->n_grp; m++) {
1570  if (tmp_grp[m].num_of_item > 0) free(tmp_grp[m].item);
1571  }
1572  free(tmp_grp);
1573  if (mynode == 0) fprintf(stderr, "Finish generating new surf_grp inf.\n");
1574  }
1575 
1576  if (tmp_int_send != NULL) free(tmp_int_send);
1577  if (tmp2_int_send != NULL) free(tmp2_int_send);
1578 
1579  /* recv_index_int=0;
1580  for(i=0;i<recv_elem_num[pesize];i++)
1581  recv_index_int+=tmp_int_recv[i];
1582  recv2_index_int=0;
1583  for(i=0;i<recv_parent_num[pesize];i++)
1584  recv2_index_int+=tmp2_int_recv[i];
1585 */
1586  if (tmp_int_recv != NULL) free(tmp_int_recv);
1587  if (tmp2_int_recv != NULL) free(tmp2_int_recv);
1588  tmp_int_send = (int *)calloc(send_ptr_num[pesize], sizeof(int));
1589  tmp_int_nodeid = (int *)calloc(recv_ptr_num[pesize], sizeof(int));
1590 
1591  if ((tmp_int_send == NULL) || (tmp_int_nodeid == NULL))
1592  HECMW_dlb_memory_exit("tmp_int_send, tmp_int_nodeid for ptr_elem sending");
1593  /* for(i=0;i<pesize+1;i++)
1594  fprintf(stderr, "vtxdist=== %d ", vtxdist[i]);
1595  fprintf(stderr, "\n");
1596  */
1597  tmp_int = 0;
1598  for (i = 0; i < send_elem_num[pesize]; i++) {
1599  for (j = mesh->elem_node_index[send_elem[i]];
1600  j < mesh->elem_node_index[send_elem[i] + 1]; j++) {
1601  tmp_int_send[tmp_int] =
1602  mesh->node_ID[(mesh->elem_node_item[j] - 1) * 2] - 1 +
1603  vtxdist[mesh->node_ID[(mesh->elem_node_item[j] - 1) * 2 + 1]];
1604  tmp_int++;
1605  }
1606  }
1607 
1608  int2_whole_send_recv(send_ptr_num[pesize], recv_ptr_num[pesize], pesize,
1609  recv_ptr_num, send_ptr_num, tmp_int_send, tmp_int_nodeid,
1610  mesh->HECMW_COMM, mynode);
1611  tmp2_int_send = (int *)calloc(send_ptr_parent_num[pesize], sizeof(int));
1612  tmp2_int_nodeid = (int *)calloc(recv_ptr_parent_num[pesize], sizeof(int));
1613 
1614  if ((tmp2_int_send == NULL) || (tmp2_int_nodeid == NULL))
1615  HECMW_dlb_memory_exit("tmp_int_send, tmp_int_nodeid for ptr_elem sending");
1616  tmp_int = 0;
1617  for (i = 0; i < send_parent_num[pesize]; i++) {
1618  for (j = mesh->elem_node_index[send_parent[i]];
1619  j < mesh->elem_node_index[send_parent[i] + 1]; j++) {
1620  tmp2_int_send[tmp_int] =
1621  mesh->node_ID[(mesh->elem_node_item[j] - 1) * 2] - 1 +
1622  vtxdist[mesh->node_ID[(mesh->elem_node_item[j] - 1) * 2 + 1]];
1623  tmp_int++;
1624  }
1625  }
1626 
1627  int2_whole_send_recv(send_ptr_parent_num[pesize], recv_ptr_parent_num[pesize],
1628  pesize, recv_ptr_parent_num, send_ptr_parent_num,
1629  tmp2_int_send, tmp2_int_nodeid, mesh->HECMW_COMM,
1630  mynode);
1631  free(tmp_int_send);
1632  free(tmp2_int_send);
1633 
1635 
1636  /* --------- start to find import_nodes according to tmp_int_nodeid,
1637  * tmp_int_peid, tmp_int_repart -------*/
1638  /* first building global_index according to node redistribution */
1639  send_node_num = (int *)calloc(pesize + 1, sizeof(int));
1640  count_node = (int *)calloc(pesize, sizeof(int));
1641  send_node = (int *)calloc(mesh->nn_internal, sizeof(int));
1642  recv_node_num = (int *)calloc(pesize + 1, sizeof(int));
1643  if ((send_node_num == NULL) || (count_node == NULL) || (send_node == NULL) ||
1644  (recv_node_num == NULL))
1646  "send_node_num, count_node, send_node, and recv_node_num");
1647  for (i = 0; i < pesize + 1; i++) {
1648  send_node_num[i] = 0;
1649  }
1650  for (i = 0; i < mesh->nn_internal; i++) {
1651  send_node_num[result->part[i] + 1]++;
1652  }
1653  for (i = 1; i < pesize + 1; i++)
1654  send_node_num[i] = send_node_num[i - 1] + send_node_num[i];
1655  for (i = 0; i < pesize; i++) count_node[i] = 0;
1656  for (i = 0; i < mesh->nn_internal; i++) {
1657  send_node[send_node_num[result->part[i]] + count_node[result->part[i]]] = i;
1658  count_node[result->part[i]]++;
1659  }
1660  tmp_node = (int *)calloc(send_node_num[pesize] + 1, sizeof(int));
1661  for (i = 0; i < send_node_num[pesize]; i++)
1662  tmp_node[i] = vtxdist[mynode] + send_node[i];
1663 
1664  stack_whole_send_recv(pesize, send_node_num, recv_node_num, mesh->HECMW_COMM,
1665  mynode);
1666  recv_node = (int *)calloc(recv_node_num[pesize] + 1, sizeof(int));
1667  if (recv_node == NULL) HECMW_dlb_memory_exit("recv_elem");
1668  int2_whole_send_recv(send_node_num[pesize], recv_node_num[pesize], pesize,
1669  recv_node_num, send_node_num, tmp_node, recv_node,
1670  mesh->HECMW_COMM, mynode);
1671  global_index = (int *)calloc(result->t_node, sizeof(int));
1672  global_index_hit = (int *)calloc(result->t_node, sizeof(int));
1673  if ((global_index == NULL) || (global_index_hit == NULL))
1674  HECMW_dlb_memory_exit("global_index");
1675  for (i = 0; i < result->t_node; i++) {
1676  global_index[i] = -1;
1677  global_index_hit[i] = -1;
1678  }
1679  for (i = 0; i < recv_node_num[pesize]; i++) {
1680  if (recv_node[i] >= result->t_node)
1681  HECMW_dlb_print_exit("!!!! wrong in recv_node");
1682  global_index[recv_node[i]] = mynode * result->t_node + i;
1683  global_index_hit[recv_node[i]] = i;
1684  }
1685  free(tmp_node);
1686 
1687  /*
1688  if(mynode>=0) {
1689  */
1690  if (mynode == 0) {
1691  tmp_recv = (int *)calloc(result->t_node, sizeof(int));
1692  if (tmp_recv == NULL) HECMW_dlb_memory_exit("tmp_recv");
1693  for (i = 1; i < pesize; i++) {
1694  HECMW_Recv(tmp_recv, result->t_node, HECMW_INT, i, HECMW_ANY_TAG,
1695  mesh->HECMW_COMM, &stat);
1696  for (j = 0; j < result->t_node; j++) {
1697  if (tmp_recv[j] >= 0) global_index[j] = tmp_recv[j];
1698  }
1699  }
1700  free(tmp_recv);
1701  for (i = 1; i < pesize; i++)
1702  HECMW_Send(global_index, result->t_node, HECMW_INT, i, 0,
1703  mesh->HECMW_COMM);
1704  } else {
1705  HECMW_Send(global_index, result->t_node, HECMW_INT, 0, 0, mesh->HECMW_COMM);
1706  HECMW_Recv(global_index, result->t_node, HECMW_INT, 0, HECMW_ANY_TAG,
1707  mesh->HECMW_COMM, &stat);
1708  }
1709 
1710  /* for(i=0;i<result->t_node;i++)
1711  fprintf(test_fp, "%d\n", global_index[i]);
1712  fclose(test_fp);
1713  */
1714  import_index = (int *)calloc(pesize + 1, sizeof(int));
1715  if (import_index == NULL) HECMW_dlb_memory_exit("import_index");
1716  for (i = 0; i < pesize + 1; i++) import_index[i] = 0;
1717  for (i = 0; i < recv_ptr_num[pesize]; i++) {
1718  tmp_int = tmp_int_nodeid[i];
1719  if (tmp_int >= result->t_node) {
1720  fprintf(stderr, "There is something wrong with data: i=%d tmp_int=%d\n",
1721  i, tmp_int);
1722  HECMW_dlb_print_exit("Please check again");
1723  }
1724  if (global_index_hit[tmp_int] == -1) {
1725  global_index_hit[tmp_int] = -2;
1726  m = global_index[tmp_int] / result->t_node;
1727  import_index[m + 1]++;
1728  }
1729  }
1730  for (i = 0; i < recv_ptr_parent_num[pesize]; i++) {
1731  tmp_int = tmp2_int_nodeid[i];
1732  if (tmp_int >= result->t_node) {
1733  fprintf(stderr, "There is something wrong with data: i=%d tmp_int=%d\n",
1734  i, tmp_int);
1735  HECMW_dlb_print_exit("Please check again");
1736  }
1737  if (global_index_hit[tmp_int] == -1) {
1738  global_index_hit[tmp_int] = -2;
1739  m = global_index[tmp_int] / result->t_node;
1740  import_index[m + 1]++;
1741  }
1742  }
1743  for (i = 1; i < pesize + 1; i++)
1744  import_index[i] = import_index[i - 1] + import_index[i];
1745  count_index = (int *)calloc(pesize, sizeof(int));
1746  for (i = 0; i < pesize; i++) count_index[i] = 0;
1747  for (i = 0; i < result->t_node; i++) {
1748  global_index_hit[i] = -1;
1749  }
1750  for (i = 0; i < recv_node_num[pesize]; i++) {
1751  global_index_hit[recv_node[i]] = i;
1752  }
1753  free(recv_node);
1754  /* ********original mesh free
1755  free(mesh->node_id);
1756 
1757 
1758  */
1759  new_mesh->n_node = recv_node_num[pesize] + import_index[pesize];
1761  (int *)calloc(new_mesh->n_dof_grp + 1, sizeof(int));
1762  if (new_mesh->node_dof_index == NULL)
1763  HECMW_dlb_memory_exit("new_mesh: node_dof_index");
1764  new_mesh->node_dof_index[0] = 0;
1766  new_mesh->node_dof_item = (int *)calloc(new_mesh->n_dof_grp, sizeof(int));
1767  if (new_mesh->node_dof_item == NULL)
1768  HECMW_dlb_memory_exit("new_mesh: node_dof_item");
1769  for (i = 0; i < new_mesh->n_dof_grp; i++)
1771 
1772  new_mesh->nn_internal = recv_node_num[pesize];
1773  new_mesh->n_dof = mesh->n_dof;
1775  new_mesh->node_ID = (int *)calloc(2 * new_mesh->n_node, sizeof(int));
1776  if (new_mesh->node_ID == NULL) HECMW_dlb_memory_exit("new_mesh: node_ID");
1777  for (i = 0; i < recv_node_num[pesize]; i++) {
1778  new_mesh->node_ID[i * 2 + 1] = mynode;
1779  new_mesh->node_ID[i * 2] = i + 1;
1780  }
1781  if ((recv_ptr_num[pesize] + recv_ptr_parent_num[pesize]) > 0) {
1782  new_mesh->elem_node_item = (int *)calloc(
1783  recv_ptr_num[pesize] + recv_ptr_parent_num[pesize], sizeof(int));
1784 
1785  if (new_mesh->elem_node_item == NULL)
1786  HECMW_dlb_memory_exit("new_mesh: elem_node_item");
1787  for (i = 0; i < recv_ptr_num[pesize]; i++) {
1788  tmp_int = tmp_int_nodeid[i];
1789  if (global_index_hit[tmp_int] != -1) {
1790  new_mesh->elem_node_item[i] = global_index_hit[tmp_int];
1791  } else {
1792  m = global_index[tmp_int] / result->t_node;
1793  global_index_hit[tmp_int] =
1794  new_mesh->nn_internal + import_index[m] + count_index[m];
1795  new_mesh->elem_node_item[i] = global_index_hit[tmp_int];
1796  new_mesh->node_ID[global_index_hit[tmp_int] * 2] =
1797  (global_index[tmp_int] % result->t_node) + 1;
1798  new_mesh->node_ID[global_index_hit[tmp_int] * 2 + 1] =
1799  global_index[tmp_int] / result->t_node;
1800  count_index[m]++;
1801  }
1802  }
1803  for (i = 0; i < recv_ptr_parent_num[pesize]; i++) {
1804  tmp_int = tmp2_int_nodeid[i];
1805  if (global_index_hit[tmp_int] != -1) {
1806  new_mesh->elem_node_item[i + recv_ptr_num[pesize]] =
1807  global_index_hit[tmp_int];
1808  } else {
1809  m = global_index[tmp_int] / result->t_node;
1810  global_index_hit[tmp_int] =
1811  new_mesh->nn_internal + import_index[m] + count_index[m];
1812  new_mesh->elem_node_item[i + recv_ptr_num[pesize]] =
1813  global_index_hit[tmp_int];
1814  new_mesh->node_ID[global_index_hit[tmp_int] * 2] =
1815  (global_index[tmp_int] % result->t_node) + 1;
1816  new_mesh->node_ID[global_index_hit[tmp_int] * 2 + 1] =
1817  global_index[tmp_int] / result->t_node;
1818  count_index[m]++;
1819  }
1820  }
1821  }
1822  for (i = 0; i < recv_ptr_num[pesize] + recv_ptr_parent_num[pesize]; i++)
1823  new_mesh->elem_node_item[i] += 1;
1824 
1826  (int *)calloc(new_mesh->ne_internal, sizeof(int));
1828  HECMW_dlb_memory_exit("new_mesh: elem_internal_list");
1829  for (i = 0; i < new_mesh->ne_internal; i++)
1830  new_mesh->elem_internal_list[i] = 0;
1831  new_mesh->ne_internal = 0;
1832  /* for(i=0;i<recv_elem_num[pesize];i++) {
1833  min_pe=mynode;
1834  for(j=new_smesh->index_elem[i];j<new_smesh->index_elem[i+1];j++)
1835  {
1836  local_nid=new_smesh->ptr_elem[j]-1;
1837  tmp_pe=new_smesh->node_id[new_smesh->n_node+local_nid];
1838  if(tmp_pe<min_pe)
1839  min_pe=tmp_pe;
1840  }
1841  if(min_pe==mynode) {
1842  new_smesh->ne_internal++;
1843  new_smesh->ne_internal_list[i]=1;
1844  }
1845  }
1846  */
1847 
1848  for (i = 0; i < new_mesh->n_elem; i++) {
1849  if (new_mesh->elem_ID[i * 2 + 1] == mynode) {
1851  new_mesh->ne_internal++;
1852  }
1853  }
1854  /*
1855  fprintf(stderr, "In pe %d ne_internal=%d n_elem=%d t_elem=%d\n", mynode,
1856  new_mesh->ne_internal, new_mesh->n_elem,
1857  t_elem);
1858  */
1859  tmp_int = 0;
1860  for (i = 1; i < pesize + 1; i++) {
1861  if ((import_index[i] - import_index[i - 1]) > 0) tmp_int++;
1862  }
1863  import_n_neighbor_pe = tmp_int;
1864  import_neighbor_pe = (int *)calloc(tmp_int, sizeof(int));
1865 
1866  tmp_int = -1;
1867  for (i = 1; i < pesize + 1; i++) {
1868  if ((import_index[i] - import_index[i - 1]) > 0) {
1869  tmp_int++;
1870  import_neighbor_pe[tmp_int] = i - 1;
1871  }
1872  }
1873  export_index = (int *)calloc(pesize + 1, sizeof(int));
1874  stack_whole_send_recv(pesize, import_index, export_index, mesh->HECMW_COMM,
1875  mynode);
1876  export_node = (int *)calloc(export_index[pesize] + 1, sizeof(int));
1877  tmp_send = (int *)calloc(import_index[pesize] + 1, sizeof(int));
1878  if (tmp_send == NULL) HECMW_dlb_memory_exit("tmp_send");
1879  for (i = 0; i < import_index[pesize]; i++)
1880  tmp_send[i] = new_mesh->node_ID[(new_mesh->nn_internal + i) * 2];
1881  /*
1882  for(i=0;i<vis_pesize+1;i++)
1883  fprintf(test_fp, "%d ", import_index[i]);
1884  fprintf(test_fp, "\n");
1885  for(i=0;i<import_index[vis_pesize];i++)
1886  fprintf(test_fp, "%d %d\n", tmp_send[i],
1887  new_smesh->node_ID[(new_smesh->nn_internal+i)*2+1]);
1888  fclose(test_fp);
1889  */
1890 
1891  int2_whole_send_recv(import_index[pesize], export_index[pesize], pesize,
1892  export_index, import_index, tmp_send, export_node,
1893  mesh->HECMW_COMM, mynode);
1894 
1895  /* for(i=0;i<pesize+1;i++)
1896  fprintf(test_fp2, "%d ", export_index[i]);
1897  fprintf(test_fp2, "\n");
1898  for(i=0;i<export_index[pesize];i++)
1899  fprintf(test_fp2, "%d \n", export_node[i]);
1900  */
1901 
1902  tmp_int = 0;
1903  for (i = 1; i < pesize + 1; i++) {
1904  if ((export_index[i] - export_index[i - 1]) > 0) tmp_int++;
1905  }
1906  export_n_neighbor_pe = tmp_int;
1907  export_neighbor_pe = (int *)calloc(tmp_int, sizeof(int));
1908  tmp_int = -1;
1909  for (i = 1; i < pesize + 1; i++) {
1910  if ((export_index[i] - export_index[i - 1]) > 0) {
1911  tmp_int++;
1912  export_neighbor_pe[tmp_int] = i - 1;
1913  }
1914  }
1915 
1916  if (export_n_neighbor_pe > import_n_neighbor_pe) {
1917  new_mesh->n_neighbor_pe = export_n_neighbor_pe;
1918  new_mesh->neighbor_pe = (int *)calloc(new_mesh->n_neighbor_pe, sizeof(int));
1919  for (i = 0; i < new_mesh->n_neighbor_pe; i++)
1920  new_mesh->neighbor_pe[i] = export_neighbor_pe[i];
1921  } else {
1922  new_mesh->n_neighbor_pe = import_n_neighbor_pe;
1923  new_mesh->neighbor_pe = (int *)calloc(new_mesh->n_neighbor_pe, sizeof(int));
1924  for (i = 0; i < new_mesh->n_neighbor_pe; i++)
1925  new_mesh->neighbor_pe[i] = import_neighbor_pe[i];
1926  }
1927  /* free(import_neighbor_pe);
1928  free(export_neighbor_pe);
1929 */
1931  (int *)calloc(new_mesh->n_neighbor_pe + 1, sizeof(int));
1933  (int *)calloc(new_mesh->n_neighbor_pe + 1, sizeof(int));
1934  new_mesh->export_index[0] = 0;
1935  new_mesh->import_index[0] = 0;
1936  for (i = 0; i < new_mesh->n_neighbor_pe; i++) {
1937  new_mesh->export_index[i + 1] = export_index[new_mesh->neighbor_pe[i] + 1];
1938  new_mesh->import_index[i + 1] = import_index[new_mesh->neighbor_pe[i] + 1];
1939  }
1940  if (import_index[pesize] > 0) {
1941  new_mesh->import_item = (int *)calloc(import_index[pesize], sizeof(int));
1942  if (new_mesh->import_item == NULL)
1943  HECMW_dlb_memory_exit("new_mesh: import_item");
1944  }
1945  for (i = 0; i < import_index[pesize]; i++)
1946  new_mesh->import_item[i] = new_mesh->nn_internal + i + 1;
1947  if (export_index[pesize] > 0) {
1948  new_mesh->export_item = (int *)calloc(export_index[pesize], sizeof(int));
1949  if (new_mesh->export_item == NULL)
1950  HECMW_dlb_memory_exit("new_mesh: export_item");
1951  }
1952  for (i = 0; i < export_index[pesize]; i++)
1953  new_mesh->export_item[i] = export_node[i];
1954  free(global_index_hit);
1955 
1957 
1958  /*
1959  }
1960  */ /* end of if(mynode in VIS_COMM) */
1961 
1962  /* ----------------sending node information --------------------- */
1963 
1964  if (new_mesh->n_node > 0) {
1965  new_mesh->node = (double *)calloc(3 * new_mesh->n_node, sizeof(double));
1966  if (new_mesh->node == NULL) HECMW_dlb_memory_exit("new_mesh: node");
1967  }
1968 
1969  tmp_node_d = (double *)calloc(send_node_num[pesize] + 1, sizeof(double));
1970  recv_node_d = (double *)calloc(recv_node_num[pesize] + 1, sizeof(double));
1971  tmp2_node_d = (double *)calloc(export_index[pesize] + 1, sizeof(double));
1972  recv2_node_d = (double *)calloc(import_index[pesize] + 1, sizeof(double));
1973  if ((tmp_node_d == NULL) || (recv_node_d == NULL) || (tmp2_node_d == NULL) ||
1974  (recv2_node_d == NULL))
1975  HECMW_dlb_memory_exit("new_mesh: recv_node");
1976  for (j = 0; j < 3; j++) {
1977  for (i = 0; i < recv_node_num[pesize]; i++) recv_node_d[i] = 0.0;
1978  for (i = 0; i < send_node_num[pesize]; i++)
1979  tmp_node_d[i] = mesh->node[send_node[i] * 3 + j];
1980  double2_whole_send_recv(send_node_num[pesize], recv_node_num[pesize],
1981  pesize, recv_node_num, send_node_num, tmp_node_d,
1982  recv_node_d, mesh->HECMW_COMM, mynode);
1983  for (i = 0; i < export_index[pesize]; i++)
1984  tmp2_node_d[i] = recv_node_d[export_node[i] - 1];
1985 
1986  double2_whole_send_recv(export_index[pesize], import_index[pesize], pesize,
1987  import_index, export_index, tmp2_node_d,
1988  recv2_node_d, mesh->HECMW_COMM, mynode);
1989  for (i = 0; i < new_mesh->nn_internal; i++)
1990  new_mesh->node[i * 3 + j] = recv_node_d[i];
1991  for (i = 0; i < import_index[pesize]; i++)
1992  new_mesh->node[(i + new_mesh->nn_internal) * 3 + j] = recv2_node_d[i];
1993  }
1994 
1995  /*
1996  fprintf(stderr, "n_node=%d nn_internal=%d recv_node_num=%d
1997  import_index_num=%d\n", new_mesh->n_node, new_mesh->nn_internal,
1998  recv_node_num[pesize],
1999  new_mesh->import_index[new_mesh->n_neighbor_pe]);
2000  */
2001  /*
2002  global_comm_table->send_node_num=send_node_num;
2003  global_comm_table->recv_node_num=recv_node_num;
2004 
2005  global_comm_table->send_node=send_node;
2006  if(mynode>=mesh_pesize) {
2007  global_comm_table->import_index=import_index;
2008  global_comm_table->export_index=export_index;
2009  global_comm_table->export_node=export_node;
2010  }
2011  */
2012  /* free(node->data);
2013  */
2014  free(tmp_int_nodeid);
2015  free(tmp2_int_nodeid);
2016 
2017  /* --------------------start sending result data ---------------------- */
2020  if (new_data->nn_component > 0) {
2021  new_data->nn_dof = (int *)calloc(new_data->nn_component, sizeof(int));
2022  new_data->node_label =
2023  (char **)calloc(new_data->nn_component, sizeof(char *));
2024  for (i = 0; i < new_data->nn_component; i++)
2025  new_data->node_label[i] = (char *)calloc(128, sizeof(char));
2026  if (new_data->nn_dof == NULL) HECMW_dlb_memory_exit("new_data: nn_dof");
2027  for (i = 0; i < new_data->nn_component; i++) {
2028  new_data->nn_dof[i] = data->nn_dof[i];
2029  strcpy(new_data->node_label[i], data->node_label[i]);
2030  }
2031  }
2032  if (new_data->ne_component > 0) {
2033  new_data->ne_dof = (int *)calloc(new_data->ne_component, sizeof(int));
2034  new_data->elem_label =
2035  (char **)calloc(new_data->ne_component, sizeof(char *));
2036  for (i = 0; i < new_data->ne_component; i++)
2037  new_data->elem_label[i] = (char *)calloc(128, sizeof(char));
2038  if (new_data->ne_dof == NULL) HECMW_dlb_memory_exit("new_data: ne_dof");
2039  for (i = 0; i < new_data->ne_component; i++) {
2040  new_data->ne_dof[i] = data->ne_dof[i];
2041  strcpy(new_data->elem_label[i], data->elem_label[i]);
2042  }
2043  }
2044  if (new_data->nn_component > 0) {
2045  tn_component = 0;
2046  for (i = 0; i < new_data->nn_component; i++)
2047  tn_component += new_data->nn_dof[i];
2049  (double *)calloc(tn_component * new_mesh->n_node, sizeof(double));
2050  if (new_data->node_val_item == NULL)
2051  HECMW_dlb_memory_exit("new_data: node_val_item");
2052  for (j = 0; j < tn_component; j++) {
2053  for (i = 0; i < send_node_num[pesize]; i++)
2054  tmp_node_d[i] = data->node_val_item[send_node[i] * tn_component + j];
2055 
2056  double2_whole_send_recv(send_node_num[pesize], recv_node_num[pesize],
2057  pesize, recv_node_num, send_node_num, tmp_node_d,
2058  recv_node_d, mesh->HECMW_COMM, mynode);
2059  for (i = 0; i < export_index[pesize]; i++)
2060  tmp2_node_d[i] = recv_node_d[export_node[i] - 1];
2061 
2062  double2_whole_send_recv(export_index[pesize], import_index[pesize],
2063  pesize, import_index, export_index, tmp2_node_d,
2064  recv2_node_d, mesh->HECMW_COMM, mynode);
2065  for (i = 0; i < new_mesh->nn_internal; i++)
2066  new_data->node_val_item[i * tn_component + j] = recv_node_d[i];
2067  for (i = 0; i < import_index[pesize]; i++)
2068  new_data
2069  ->node_val_item[(i + new_mesh->nn_internal) * tn_component + j] =
2070  recv2_node_d[i];
2071  }
2072  }
2073 
2074  free(tmp_node_d);
2075  free(recv_node_d);
2076  free(tmp2_node_d);
2077  free(recv2_node_d);
2079  (int *)calloc(new_mesh->n_node, sizeof(int));
2081  HECMW_dlb_memory_exit("new_mesh: when_i_was_refined_node");
2082  tmp_node_i = (int *)calloc(send_node_num[pesize] + 1, sizeof(int));
2083  recv_node_i = (int *)calloc(recv_node_num[pesize] + 1, sizeof(int));
2084  tmp2_node_i = (int *)calloc(export_index[pesize] + 1, sizeof(int));
2085  recv2_node_i = (int *)calloc(import_index[pesize] + 1, sizeof(int));
2086  if ((tmp_node_i == NULL) || (recv_node_i == NULL))
2087  HECMW_dlb_memory_exit("new_mesh: when_i_was_refined_node");
2088  for (i = 0; i < recv_node_num[pesize]; i++) recv_node_i[i] = -1;
2089  for (i = 0; i < send_node_num[pesize]; i++)
2090  tmp_node_i[i] = mesh->when_i_was_refined_node[send_node[i]];
2091 
2092  int2_whole_send_recv(send_node_num[pesize], recv_node_num[pesize], pesize,
2093  recv_node_num, send_node_num, tmp_node_i, recv_node_i,
2094  mesh->HECMW_COMM, mynode);
2095  for (i = 0; i < export_index[pesize]; i++)
2096  tmp2_node_i[i] = recv_node_i[export_node[i] - 1];
2097 
2098  int2_whole_send_recv(export_index[pesize], import_index[pesize], pesize,
2099  import_index, export_index, tmp2_node_i, recv2_node_i,
2100  mesh->HECMW_COMM, mynode);
2101  for (i = 0; i < new_mesh->nn_internal; i++)
2102  new_mesh->when_i_was_refined_node[i] = recv_node_i[i];
2103  for (i = 0; i < import_index[pesize]; i++)
2105  recv2_node_i[i];
2106 
2107  if (new_mesh->node_group->n_grp > 0) {
2108  tmp_grp =
2109  (Tmp_grp_inf *)calloc(new_mesh->node_group->n_grp, sizeof(Tmp_grp_inf));
2110  if (tmp_grp == NULL) HECMW_dlb_memory_exit("tmp_grp");
2111 
2112  for (m = 0; m < new_mesh->node_group->n_grp; m++)
2113  tmp_grp[m].num_of_item = 0;
2114  for (m = 0; m < new_mesh->node_group->n_grp; m++) {
2115  if (m == 0) {
2116  tmp_elem_grp = (int *)calloc(mesh->nn_internal, sizeof(int));
2117  if (tmp_elem_grp == NULL) HECMW_dlb_memory_exit("tmp_elem_grp");
2118  }
2119  for (i = 0; i < mesh->nn_internal; i++) tmp_elem_grp[i] = 0;
2120  if ((mesh->node_group->grp_index[m + 1] -
2121  mesh->node_group->grp_index[m]) > 0) {
2122  for (i = mesh->node_group->grp_index[m];
2123  i < mesh->node_group->grp_index[m + 1]; i++) {
2124  if (mesh->node_group->grp_item[i] <= mesh->nn_internal)
2125  tmp_elem_grp[mesh->node_group->grp_item[i] - 1] = 1;
2126  }
2127  }
2128  for (i = 0; i < recv_node_num[pesize]; i++) recv_node_i[i] = -1;
2129  for (i = 0; i < send_node_num[pesize]; i++)
2130  tmp_node_i[i] = tmp_elem_grp[send_node[i]];
2131 
2132  int2_whole_send_recv(send_node_num[pesize], recv_node_num[pesize], pesize,
2133  recv_node_num, send_node_num, tmp_node_i,
2134  recv_node_i, mesh->HECMW_COMM, mynode);
2135  for (i = 0; i < export_index[pesize]; i++)
2136  tmp2_node_i[i] = recv_node_i[export_node[i] - 1];
2137 
2138  int2_whole_send_recv(export_index[pesize], import_index[pesize], pesize,
2139  import_index, export_index, tmp2_node_i,
2140  recv2_node_i, mesh->HECMW_COMM, mynode);
2141  num_grp_item = 0;
2142  for (i = 0; i < new_mesh->nn_internal; i++) {
2143  if (recv_node_i[i] == 1) num_grp_item++;
2144  }
2145  for (i = 0; i < import_index[pesize]; i++) {
2146  if (recv2_node_i[i] == 1) num_grp_item++;
2147  }
2148  tmp_grp[m].num_of_item = num_grp_item;
2149  if (num_grp_item > 0) {
2150  tmp_grp[m].item = (int *)calloc(num_grp_item, sizeof(int));
2151  if (tmp_grp[m].item == NULL) HECMW_dlb_memory_exit("tmp_grp:item");
2152  tmp_int = 0;
2153  for (i = 0; i < new_mesh->nn_internal; i++) {
2154  if (recv_node_i[i] == 1) {
2155  tmp_grp[m].item[tmp_int] = i + 1;
2156  tmp_int++;
2157  }
2158  }
2159  for (i = 0; i < import_index[pesize]; i++) {
2160  if (recv2_node_i[i] == 1) {
2161  tmp_grp[m].item[tmp_int] = new_mesh->nn_internal + i + 1;
2162  tmp_int++;
2163  }
2164  }
2165  }
2166  }
2167 
2168  free(tmp_elem_grp);
2169  num_grp_item = 0;
2170  for (m = 0; m < new_mesh->node_group->n_grp; m++)
2171  num_grp_item += tmp_grp[m].num_of_item;
2173  (int *)calloc(new_mesh->node_group->n_grp + 1, sizeof(int));
2174  if (new_mesh->node_group->grp_index == NULL)
2175  HECMW_dlb_memory_exit("new_mesh: node_grp: grp_index");
2176  new_mesh->node_group->grp_index[0] = 0;
2177  for (m = 0; m < new_mesh->node_group->n_grp; m++)
2178  new_mesh->node_group->grp_index[m + 1] =
2179  new_mesh->node_group->grp_index[m] + tmp_grp[m].num_of_item;
2180  if (num_grp_item > 0) {
2181  new_mesh->node_group->grp_item = (int *)calloc(num_grp_item, sizeof(int));
2182  if (new_mesh->node_group->grp_item == NULL)
2183  HECMW_dlb_memory_exit("new_mesh: node_grp: grp_item");
2184  tmp_int = 0;
2185  for (m = 0; m < new_mesh->node_group->n_grp; m++) {
2186  for (i = 0; i < tmp_grp[m].num_of_item; i++) {
2187  new_mesh->node_group->grp_item[tmp_int] = tmp_grp[m].item[i];
2188  tmp_int++;
2189  }
2190  }
2191  }
2192  for (m = 0; m < new_mesh->node_group->n_grp; m++) {
2193  if (tmp_grp[m].num_of_item > 0) free(tmp_grp[m].item);
2194  }
2195  free(tmp_grp);
2196  if (mynode == 0) fprintf(stderr, "Finish generating new node_grp inf.\n");
2197  }
2198 
2199  free(tmp_node_i);
2200  free(tmp2_node_i);
2201  free(recv_node_i);
2202  free(recv2_node_i);
2203 
2204  nvtxs = new_mesh->nn_internal;
2205  new_vtxdist = (int *)calloc(pesize + 1, sizeof(int));
2206 
2207  if (mynode == 0) {
2208  new_vtxdist[0] = 0;
2209  new_vtxdist[1] = nvtxs;
2210  tmp_sum = nvtxs;
2211  for (i = 1; i < pesize; i++) {
2212  HECMW_Recv(&tmp_nvtxs, 1, HECMW_INT, i, HECMW_ANY_TAG, mesh->HECMW_COMM,
2213  &stat);
2214  tmp_sum += tmp_nvtxs;
2215  new_vtxdist[i + 1] = tmp_sum;
2216  }
2217  for (i = 1; i < pesize; i++)
2218  HECMW_Send(new_vtxdist, pesize + 1, HECMW_INT, i, 0, mesh->HECMW_COMM);
2219  } else {
2220  HECMW_Send(&nvtxs, 1, HECMW_INT, 0, 0, mesh->HECMW_COMM);
2221  HECMW_Recv(new_vtxdist, pesize + 1, HECMW_INT, 0, HECMW_ANY_TAG,
2222  mesh->HECMW_COMM, &stat);
2223  }
2224  new_mesh->global_node_ID = (int *)calloc(new_mesh->n_node, sizeof(int));
2225  if (new_mesh->global_node_ID == NULL)
2226  HECMW_dlb_print_exit("new_mesh: global_node_ID");
2227  for (i = 0; i < new_mesh->n_node; i++)
2228  new_mesh->global_node_ID[i] =
2229  new_vtxdist[new_mesh->node_ID[i * 2 + 1]] + new_mesh->node_ID[i * 2];
2230  free(new_vtxdist);
2231  /*
2232  for(i=0;i<result->t_node;i++)
2233  global_new2old[i]=-1;
2234 
2235  for(i=0; i<result->t_node;i++) {
2236  tmp_peid=global_index[i] / result->t_node;
2237  tmp_lid=global_index[i] % result->t_node;
2238  global_new2old[new_vtxdist[tmp_peid]+tmp_lid]=i;
2239  }
2240  */
2241  free(global_index);
2242  /* set new data structure */
2243  /* if(mynode>=mesh_pesize) {
2244  new_data->node_val_item=(double *)calloc(new_smesh->n_node,
2245  sizeof(double));
2246  if(new_data->node_val_item==NULL)
2247  HECMW_dlb_memory_exit("new_data: node_val_item");
2248  }
2249  */
2250 
2251  if (new_mesh->n_elem_type > 1) {
2252  new2old = (int *)calloc(new_mesh->n_elem, sizeof(int));
2253  if (new2old == NULL) HECMW_dlb_memory_exit("new2old");
2255  (int *)calloc(new_mesh->n_elem_type + 1, sizeof(int));
2256  if (new_mesh->elem_type_index == NULL)
2257  HECMW_dlb_memory_exit("new_mesh: elem_type_index");
2258  for (i = 0; i < new_mesh->n_elem_type + 1; i++)
2259  new_mesh->elem_type_index[i] = 0;
2261  (int *)calloc(new_mesh->n_elem_type, sizeof(int));
2262  if (new_mesh->elem_type_item == NULL)
2263  HECMW_dlb_memory_exit("new_mesh: elem_type_item");
2264  for (i = 0; i < new_mesh->n_elem_type; i++)
2266 
2267  for (i = 0; i < new_mesh->n_elem; i++) {
2268  for (j = 0; j < new_mesh->n_elem_type; j++) {
2269  if (new_mesh->elem_type[i] == new_mesh->elem_type_item[j]) {
2270  new_mesh->elem_type_index[j + 1]++;
2271  break;
2272  }
2273  }
2274  }
2275  for (j = 1; j < new_mesh->n_elem_type + 1; j++)
2277  /*
2278  fprintf(stderr, "new_mesh: elem_type_index= %d %d %d\n",
2279  new_mesh->elem_type_index[0], new_mesh->elem_type_index[1],
2280  new_mesh->elem_type_index[2]);
2281  */
2282  count_elem_index = (int *)calloc(new_mesh->n_elem_type, sizeof(int));
2283  if (count_elem_index == NULL)
2284  HECMW_dlb_memory_exit("new_mesh: count_elem_index");
2285  for (i = 0; i < new_mesh->n_elem_type; i++) count_elem_index[i] = 0;
2286  for (i = 0; i < new_mesh->n_elem; i++) {
2287  for (j = 0; j < new_mesh->n_elem_type; j++) {
2288  if (new_mesh->elem_type[i] == new_mesh->elem_type_item[j]) {
2289  new2old[i] = new_mesh->elem_type_index[j] + count_elem_index[j];
2290  count_elem_index[j]++;
2291  break;
2292  }
2293  }
2294  }
2295  free(count_elem_index);
2296  new_tmp = (int *)calloc(new_mesh->n_elem, sizeof(int));
2297  if (new_tmp == NULL) HECMW_dlb_memory_exit("new_tmp");
2298  for (i = 0; i < new_mesh->n_elem; i++)
2299  new_tmp[i] = new_mesh->elem_type[new2old[i]];
2300  free(new_mesh->elem_type);
2301  new_mesh->elem_type = new_tmp;
2302  new_tmp = (int *)calloc(new_mesh->n_elem, sizeof(int));
2303  if (new_tmp == NULL) HECMW_dlb_memory_exit("new_tmp");
2304  for (i = 0; i < new_mesh->n_elem; i++)
2305  new_tmp[i] = new_mesh->section_ID[new2old[i]];
2306  free(new_mesh->section_ID);
2307  new_mesh->section_ID = new_tmp;
2308  new_tmp = (int *)calloc(new_mesh->n_elem, sizeof(int));
2309  if (new_tmp == NULL) HECMW_dlb_memory_exit("new_tmp");
2310  for (i = 0; i < new_mesh->n_elem; i++)
2311  new_tmp[i] = new_mesh->elem_mat_ID_item[new2old[i]];
2312  free(new_mesh->elem_mat_ID_item);
2313  new_mesh->elem_mat_ID_item = new_tmp;
2314 
2315  /*
2316  fprintf(test_fp, "n_node=%d n_elem=%d nn_internal=%d ne_internal=%d\n",
2317  new_mesh->n_node, new_mesh->n_elem,
2318  new_mesh->nn_internal, new_mesh->ne_internal);
2319  for(i=0;i<new_mesh->n_elem;i++)
2320  fprintf(test_fp, "%d %d\n", i,
2321  new_mesh->elem_node_index[i+1]-new_mesh->elem_node_index[i]);
2322  fclose(test_fp);
2323  */
2324 
2325  new_tmp2 = (int *)calloc(new_mesh->n_elem + 1, sizeof(int));
2326  for (i = 0; i < new_mesh->n_elem + 1; i++)
2327  new_tmp2[i] = new_mesh->elem_node_index[i];
2328  new_tmp = (int *)calloc(new_mesh->n_elem, sizeof(int));
2329  if (new_tmp == NULL) HECMW_dlb_memory_exit("new_tmp");
2330  for (i = 0; i < new_mesh->n_elem; i++)
2331  new_tmp[i] = new_mesh->elem_node_index[new2old[i] + 1] -
2332  new_mesh->elem_node_index[new2old[i]];
2333  for (i = 1; i < new_mesh->n_elem + 1; i++)
2335  new_mesh->elem_node_index[i - 1] + new_tmp[i - 1];
2336  free(new_tmp);
2337  new_tmp =
2338  (int *)calloc(new_mesh->elem_node_index[new_mesh->n_elem], sizeof(int));
2339  if (new_tmp == NULL) HECMW_dlb_memory_exit("new_tmp");
2340  for (i = 0; i < new_mesh->n_elem; i++) {
2341  for (j = new_tmp2[new2old[i]]; j < new_tmp2[new2old[i] + 1]; j++) {
2342  new_tmp[new_mesh->elem_node_index[i] + j - new_tmp2[new2old[i]]] =
2344  }
2345  }
2346  free(new_mesh->elem_node_item);
2347  new_mesh->elem_node_item = new_tmp;
2348  free(new_tmp2);
2349  new_tmp = (int *)calloc(new_mesh->n_elem * 2, sizeof(int));
2350  if (new_tmp == NULL) HECMW_dlb_memory_exit("new_tmp");
2351  for (i = 0; i < new_mesh->n_elem; i++) {
2352  new_tmp[i * 2] = new_mesh->elem_ID[new2old[i] * 2];
2353  new_tmp[i * 2 + 1] = new_mesh->elem_ID[new2old[i] * 2 + 1];
2354  }
2355  free(new_mesh->elem_ID);
2356  new_mesh->elem_ID = new_tmp;
2357  new_tmp = (int *)calloc(new_mesh->n_elem, sizeof(int));
2358  if (new_tmp == NULL) HECMW_dlb_memory_exit("new_tmp");
2359  for (i = 0; i < new_mesh->n_elem; i++)
2360  new_tmp[i] = new_mesh->global_elem_ID[new2old[i]];
2361  free(new_mesh->global_elem_ID);
2362  new_mesh->global_elem_ID = new_tmp;
2363 #ifdef TEST
2364  for (i = 0; i < new_mesh->n_elem; i++)
2365  fprintf(test_fp, "i= %d elem_ID=%d %d\n", i, new_mesh->elem_ID[i * 2],
2366  new_mesh->elem_ID[i * 2 + 1]);
2367  fclose(test_fp);
2368 #endif
2369  old2new = (int *)calloc(new_mesh->n_elem, sizeof(int));
2370  if (old2new == NULL) HECMW_dlb_memory_exit("old2new");
2371  for (i = 0; i < new_mesh->n_elem; i++) old2new[new2old[i]] = i;
2372  new_mesh->ne_internal = 0;
2373  for (i = 0; i < new_mesh->n_elem; i++) {
2374  if (new_mesh->elem_ID[i * 2 + 1] == mynode) {
2376  new_mesh->ne_internal++;
2377  }
2378  }
2379 
2381  new_mesh->zero = mesh->zero;
2382  new_mesh->PETOT = mesh->PETOT;
2384  new_mesh->errnof = mesh->errnof;
2387 
2389  (int *)calloc(new_mesh->n_neighbor_pe + 1, sizeof(int));
2390  if (new_mesh->shared_index == NULL)
2391  HECMW_dlb_memory_exit("new_mesh: shared_index");
2392  for (i = 0; i < new_mesh->n_neighbor_pe + 1; i++)
2393  new_mesh->shared_index[i] = 0;
2394 
2395  new_tmp = (int *)calloc(new_mesh->n_elem, sizeof(int));
2396  if (new_tmp == NULL) HECMW_dlb_memory_exit("new_tmp");
2397  for (i = 0; i < new_mesh->n_elem; i++)
2398  new_tmp[i] = new_mesh->when_i_was_refined_elem[new2old[i]];
2400  new_mesh->when_i_was_refined_elem = new_tmp;
2401 
2402  new_tmp = (int *)calloc(new_mesh->n_elem, sizeof(int));
2403  if (new_tmp == NULL) HECMW_dlb_memory_exit("new_tmp");
2404  for (i = 0; i < new_mesh->n_elem; i++)
2405  new_tmp[i] = new_mesh->adapt_parent_type[new2old[i]];
2406  free(new_mesh->adapt_parent_type);
2407  new_mesh->adapt_parent_type = new_tmp;
2408 
2409  new_tmp = (int *)calloc(new_mesh->n_elem, sizeof(int));
2410  if (new_tmp == NULL) HECMW_dlb_memory_exit("new_tmp");
2411  for (i = 0; i < new_mesh->n_elem; i++)
2412  new_tmp[i] = new_mesh->adapt_type[new2old[i]];
2413  free(new_mesh->adapt_type);
2414  new_mesh->adapt_type = new_tmp;
2415 
2416  new_tmp = (int *)calloc(new_mesh->n_elem, sizeof(int));
2417  if (new_tmp == NULL) HECMW_dlb_memory_exit("new_tmp");
2418  for (i = 0; i < new_mesh->n_elem; i++)
2419  new_tmp[i] = new_mesh->adapt_level[new2old[i]];
2420  free(new_mesh->adapt_level);
2421  new_mesh->adapt_level = new_tmp;
2422 
2423  new_tmp = (int *)calloc(new_mesh->n_elem * 2, sizeof(int));
2424  if (new_tmp == NULL) HECMW_dlb_memory_exit("new_tmp");
2425  for (i = 0; i < new_mesh->n_elem; i++) {
2426  new_tmp[i * 2] = new_mesh->adapt_parent[new2old[i] * 2];
2427  new_tmp[i * 2 + 1] = new_mesh->adapt_parent[new2old[i] * 2 + 1];
2428  }
2429  free(new_mesh->adapt_parent);
2430  new_mesh->adapt_parent = new_tmp;
2431 
2432  new_tmp2 = (int *)calloc(new_mesh->n_elem + 1, sizeof(int));
2433  for (i = 0; i < new_mesh->n_elem + 1; i++)
2434  new_tmp2[i] = new_mesh->adapt_children_index[i];
2435  new_tmp = (int *)calloc(new_mesh->n_elem, sizeof(int));
2436  if (new_tmp == NULL) HECMW_dlb_memory_exit("new_tmp");
2437  for (i = 0; i < new_mesh->n_elem; i++)
2438  new_tmp[i] = new_mesh->adapt_children_index[new2old[i] + 1] -
2439  new_mesh->adapt_children_index[new2old[i]];
2440  for (i = 1; i < new_mesh->n_elem + 1; i++)
2442  new_mesh->adapt_children_index[i - 1] + new_tmp[i - 1];
2443  free(new_tmp);
2444  new_tmp = (int *)calloc(
2445  new_mesh->adapt_children_index[new_mesh->n_elem] * 2, sizeof(int));
2446  if (new_tmp == NULL) HECMW_dlb_memory_exit("new_tmp");
2447  for (i = 0; i < new_mesh->n_elem; i++) {
2448  for (j = new_tmp2[new2old[i]]; j < new_tmp2[new2old[i] + 1]; j++) {
2449  new_tmp[(new_mesh->adapt_children_index[i] + j - new_tmp2[new2old[i]]) *
2450  2] = new_mesh->adapt_children_item[j * 2];
2451  new_tmp[(new_mesh->adapt_children_index[i] + j - new_tmp2[new2old[i]]) *
2452  2 +
2453  1] = new_mesh->adapt_children_item[j * 2 + 1];
2454  }
2455  }
2457  new_mesh->adapt_children_item = new_tmp;
2458  free(new_tmp2);
2459  if (mesh->section != NULL) {
2460  new_mesh->section =
2461  (struct hecmwST_section *)calloc(1, sizeof(struct hecmwST_section));
2462  if (new_mesh->section == NULL) HECMW_dlb_memory_exit("new_mesh: section");
2465  (int *)calloc(new_mesh->section->n_sect, sizeof(int));
2467  (int *)calloc(new_mesh->section->n_sect, sizeof(int));
2468  if ((new_mesh->section->sect_type == NULL) ||
2469  (new_mesh->section->sect_opt == NULL))
2470  HECMW_dlb_memory_exit("new_mesh: section");
2472  (int *)calloc(new_mesh->section->n_sect + 1, sizeof(int));
2473  new_mesh->section->sect_mat_ID_item = (int *)calloc(
2474  mesh->section->sect_mat_ID_index[mesh->section->n_sect], sizeof(int));
2475  for (i = 0; i < mesh->section->n_sect; i++)
2477  for (i = 0; i < mesh->section->n_sect; i++)
2479  for (i = 0; i < mesh->section->n_sect + 1; i++)
2482  for (i = 0; i < mesh->section->sect_mat_ID_index[mesh->section->n_sect];
2483  i++)
2487  (int *)calloc(new_mesh->section->n_sect + 1, sizeof(int));
2489  (int *)calloc(new_mesh->section->n_sect + 1, sizeof(int));
2490 
2491  for (i = 0; i < mesh->section->n_sect + 1; i++)
2493  for (i = 0; i < mesh->section->n_sect + 1; i++)
2495  new_mesh->section->sect_I_item = (int *)calloc(
2497  sizeof(int));
2498  new_mesh->section->sect_R_item = (double *)calloc(
2500  sizeof(double));
2501  for (i = 0;
2504  for (i = 0;
2507  }
2508  if (mesh->material != NULL) {
2509  new_mesh->material =
2510  (struct hecmwST_material *)calloc(1, sizeof(struct hecmwST_material));
2511  if (new_mesh->material == NULL)
2512  HECMW_dlb_memory_exit("new_mesh: material");
2518  (char **)calloc(new_mesh->material->n_mat, sizeof(char *));
2519  if (new_mesh->material->mat_name == NULL)
2520  HECMW_dlb_memory_exit("new_mesh: material");
2521  for (i = 0; i < new_mesh->material->n_mat; i++) {
2522  new_mesh->material->mat_name[i] = (char *)calloc(128, sizeof(char));
2523  sprintf(new_mesh->material->mat_name[i], "%s",
2524  mesh->material->mat_name[i]);
2525  }
2527  (int *)calloc(new_mesh->material->n_mat + 1, sizeof(int));
2529  (int *)calloc(new_mesh->material->n_mat_item + 1, sizeof(int));
2531  (int *)calloc(new_mesh->material->n_mat_subitem + 1, sizeof(int));
2532  if ((new_mesh->material->mat_item_index == NULL) ||
2535  HECMW_dlb_memory_exit("new_mesh: material");
2536  for (i = 0; i < new_mesh->material->n_mat + 1; i++)
2539  for (i = 0; i < new_mesh->material->n_mat_item + 1; i++)
2542  for (i = 0; i < new_mesh->material->n_mat_subitem + 1; i++)
2545  new_mesh->material->mat_val = (double *)calloc(
2548  sizeof(double));
2549  new_mesh->material->mat_temp = (double *)calloc(
2552  sizeof(double));
2553  if ((new_mesh->material->mat_val == NULL) ||
2554  (new_mesh->material->mat_temp == NULL))
2555  HECMW_dlb_memory_exit("new_mesh: material");
2556  for (i = 0; i < new_mesh->material
2558  i++) {
2561  }
2562  }
2563  if (mesh->mpc != NULL) {
2564  new_mesh->mpc =
2565  (struct hecmwST_mpc *)calloc(1, sizeof(struct hecmwST_mpc));
2566  if (new_mesh->mpc == NULL) HECMW_dlb_memory_exit("new_mesh: mpc");
2567  new_mesh->mpc->n_mpc = mesh->mpc->n_mpc;
2568  if (new_mesh->mpc->n_mpc != 0) {
2569  new_mesh->mpc->mpc_index =
2570  (int *)calloc(new_mesh->mpc->n_mpc + 1, sizeof(int));
2571  new_mesh->mpc->mpc_item = (int *)calloc(
2572  mesh->mpc->mpc_index[new_mesh->mpc->n_mpc], sizeof(int));
2573  new_mesh->mpc->mpc_dof = (int *)calloc(
2574  mesh->mpc->mpc_index[new_mesh->mpc->n_mpc], sizeof(int));
2575  new_mesh->mpc->mpc_val = (double *)calloc(
2576  mesh->mpc->mpc_index[new_mesh->mpc->n_mpc], sizeof(double));
2577  for (i = 0; i < mesh->mpc->n_mpc + 1; i++)
2578  new_mesh->mpc->mpc_index[i] = mesh->mpc->mpc_index[i];
2579  for (i = 0; i < mesh->mpc->n_mpc; i++) {
2580  new_mesh->mpc->mpc_item[i] = mesh->mpc->mpc_item[i];
2581  new_mesh->mpc->mpc_dof[i] = mesh->mpc->mpc_dof[i];
2582  new_mesh->mpc->mpc_val[i] = mesh->mpc->mpc_val[i];
2583  }
2584  }
2585  }
2586  if (mesh->amp != NULL) {
2587  new_mesh->amp = (struct hecmwST_amplitude *)calloc(
2588  1, sizeof(struct hecmwST_amplitude));
2589  if (new_mesh->amp == NULL) HECMW_dlb_memory_exit("new_mesh: amp");
2590  new_mesh->amp->n_amp = mesh->amp->n_amp;
2591  if (new_mesh->amp->n_amp != 0) {
2593  (int *)calloc(mesh->amp->n_amp, sizeof(int));
2595  (int *)calloc(mesh->amp->n_amp, sizeof(int));
2597  (int *)calloc(mesh->amp->n_amp, sizeof(int));
2598  new_mesh->amp->amp_index =
2599  (int *)calloc(mesh->amp->n_amp + 1, sizeof(int));
2600  new_mesh->amp->amp_val = (double *)calloc(
2601  mesh->amp->amp_index[mesh->amp->n_amp], sizeof(double));
2602  new_mesh->amp->amp_table = (double *)calloc(
2603  mesh->amp->amp_index[mesh->amp->n_amp], sizeof(double));
2604  if ((new_mesh->amp->amp_type_definition == NULL) ||
2605  (new_mesh->amp->amp_type_time == NULL) ||
2606  (new_mesh->amp->amp_type_value == NULL) ||
2607  (new_mesh->amp->amp_index == NULL) ||
2608  (new_mesh->amp->amp_val == NULL) ||
2609  (new_mesh->amp->amp_table == NULL))
2610  HECMW_dlb_memory_exit("new_mesh: amp");
2611  for (i = 0; i < mesh->amp->n_amp; i++) {
2616  }
2617  for (i = 0; i < mesh->amp->n_amp + 1; i++)
2618  new_mesh->amp->amp_index[i] = mesh->amp->amp_index[i];
2619  for (i = 0; i < mesh->amp->amp_index[mesh->amp->n_amp]; i++) {
2620  new_mesh->amp->amp_val[i] = mesh->amp->amp_val[i];
2621  new_mesh->amp->amp_table[i] = mesh->amp->amp_table[i];
2622  }
2623  }
2624  }
2625 
2626  if (new_mesh->elem_group->n_grp > 0) {
2628  new_tmp = (int *)calloc(
2630  sizeof(int));
2631  for (i = 0;
2633  i++)
2634  new_tmp[i] = old2new[new_mesh->elem_group->grp_item[i] - 1] + 1;
2635  free(new_mesh->elem_group->grp_item);
2636  new_mesh->elem_group->grp_item = new_tmp;
2637  }
2638  }
2639 
2640  if (new_mesh->surf_group->n_grp > 0) {
2642  new_tmp = (int *)calloc(
2644  sizeof(int));
2645  for (i = 0;
2647  i++) {
2648  new_tmp[i * 2] =
2649  old2new[new_mesh->surf_group->grp_item[i * 2] - 1] + 1;
2650  new_tmp[i * 2 + 1] = new_mesh->surf_group->grp_item[i * 2 + 1];
2651  }
2652  free(new_mesh->surf_group->grp_item);
2653  new_mesh->surf_group->grp_item = new_tmp;
2654  }
2655  }
2656  } else {
2657  new_mesh->elem_type_index = (int *)calloc(2, sizeof(int));
2658  if (new_mesh->elem_type_index == NULL)
2659  HECMW_dlb_memory_exit("new_mesh: elem_type_index");
2660  for (i = 0; i < new_mesh->n_elem_type + 1; i++)
2661  new_mesh->elem_type_index[i] = 0;
2663  (int *)calloc(new_mesh->n_elem_type, sizeof(int));
2664  if (new_mesh->elem_type_item == NULL)
2665  HECMW_dlb_memory_exit("new_mesh: elem_type_item");
2666  for (i = 0; i < new_mesh->n_elem_type; i++)
2669 
2670  /*
2671  fprintf(stderr, "new_mesh: elem_type_index= %d %d %d\n",
2672  new_mesh->elem_type_index[0], new_mesh->elem_type_index[1],
2673  new_mesh->elem_type_index[2]);
2674  */
2675 
2676  /*
2677  for(i=0;i<new_mesh->n_elem;i++)
2678  fprintf(test_fp, "i= %d elem_ID=%d %d\n", i,
2679  new_mesh->elem_ID[i*2], new_mesh->elem_ID[i*2+1]);
2680  fclose(test_fp);
2681  */
2683  new_mesh->zero = mesh->zero;
2684  new_mesh->PETOT = mesh->PETOT;
2686  new_mesh->errnof = mesh->errnof;
2689 
2691  (int *)calloc(new_mesh->n_neighbor_pe + 1, sizeof(int));
2692  if (new_mesh->shared_index == NULL)
2693  HECMW_dlb_memory_exit("new_mesh: shared_index");
2694  for (i = 0; i < new_mesh->n_neighbor_pe + 1; i++)
2695  new_mesh->shared_index[i] = 0;
2696 
2697  if (mesh->section != NULL) {
2698  new_mesh->section =
2699  (struct hecmwST_section *)calloc(1, sizeof(struct hecmwST_section));
2700  if (new_mesh->section == NULL) HECMW_dlb_memory_exit("new_mesh: section");
2703  (int *)calloc(new_mesh->section->n_sect, sizeof(int));
2705  (int *)calloc(new_mesh->section->n_sect, sizeof(int));
2706  if ((new_mesh->section->sect_type == NULL) ||
2707  (new_mesh->section->sect_opt == NULL))
2708  HECMW_dlb_memory_exit("new_mesh: section");
2710  (int *)calloc(new_mesh->section->n_sect + 1, sizeof(int));
2711  new_mesh->section->sect_mat_ID_item = (int *)calloc(
2712  mesh->section->sect_mat_ID_index[mesh->section->n_sect], sizeof(int));
2713  for (i = 0; i < mesh->section->n_sect; i++)
2715  for (i = 0; i < mesh->section->n_sect; i++)
2717  for (i = 0; i < mesh->section->n_sect + 1; i++)
2720  for (i = 0; i < mesh->section->sect_mat_ID_index[mesh->section->n_sect];
2721  i++)
2725  (int *)calloc(new_mesh->section->n_sect + 1, sizeof(int));
2727  (int *)calloc(new_mesh->section->n_sect + 1, sizeof(int));
2728 
2729  for (i = 0; i < mesh->section->n_sect + 1; i++)
2731  for (i = 0; i < mesh->section->n_sect + 1; i++)
2733  }
2734  if (mesh->material != NULL) {
2735  new_mesh->material =
2736  (struct hecmwST_material *)calloc(1, sizeof(struct hecmwST_material));
2737  if (new_mesh->material == NULL)
2738  HECMW_dlb_memory_exit("new_mesh: material");
2744  (char **)calloc(new_mesh->material->n_mat, sizeof(char *));
2745  if (new_mesh->material->mat_name == NULL)
2746  HECMW_dlb_memory_exit("new_mesh: material");
2747  for (i = 0; i < new_mesh->material->n_mat; i++) {
2748  new_mesh->material->mat_name[i] = (char *)calloc(128, sizeof(char));
2749  sprintf(new_mesh->material->mat_name[i], "%s",
2750  mesh->material->mat_name[i]);
2751  }
2753  (int *)calloc(new_mesh->material->n_mat + 1, sizeof(int));
2755  (int *)calloc(new_mesh->material->n_mat_item + 1, sizeof(int));
2757  (int *)calloc(new_mesh->material->n_mat_subitem + 1, sizeof(int));
2758  if ((new_mesh->material->mat_item_index == NULL) ||
2761  HECMW_dlb_memory_exit("new_mesh: material");
2762  for (i = 0; i < new_mesh->material->n_mat + 1; i++)
2765  for (i = 0; i < new_mesh->material->n_mat_item + 1; i++)
2768  for (i = 0; i < new_mesh->material->n_mat_subitem + 1; i++)
2771  new_mesh->material->mat_val = (double *)calloc(
2774  sizeof(double));
2775  new_mesh->material->mat_temp = (double *)calloc(
2778  sizeof(double));
2779  if ((new_mesh->material->mat_val == NULL) ||
2780  (new_mesh->material->mat_temp == NULL))
2781  HECMW_dlb_memory_exit("new_mesh: material");
2782  for (i = 0; i < new_mesh->material
2784  i++) {
2787  }
2788  }
2789  if (mesh->mpc != NULL) {
2790  new_mesh->mpc =
2791  (struct hecmwST_mpc *)calloc(1, sizeof(struct hecmwST_mpc));
2792  if (new_mesh->mpc == NULL) HECMW_dlb_memory_exit("new_mesh: mpc");
2793  new_mesh->mpc->n_mpc = mesh->mpc->n_mpc;
2794  if (new_mesh->mpc->n_mpc != 0) {
2795  new_mesh->mpc->mpc_index =
2796  (int *)calloc(new_mesh->mpc->n_mpc + 1, sizeof(int));
2797  new_mesh->mpc->mpc_item = (int *)calloc(
2798  mesh->mpc->mpc_index[new_mesh->mpc->n_mpc], sizeof(int));
2799  new_mesh->mpc->mpc_dof = (int *)calloc(
2800  mesh->mpc->mpc_index[new_mesh->mpc->n_mpc], sizeof(int));
2801  new_mesh->mpc->mpc_val = (double *)calloc(
2802  mesh->mpc->mpc_index[new_mesh->mpc->n_mpc], sizeof(double));
2803  for (i = 0; i < mesh->mpc->n_mpc + 1; i++)
2804  new_mesh->mpc->mpc_index[i] = mesh->mpc->mpc_index[i];
2805  for (i = 0; i < mesh->mpc->n_mpc; i++) {
2806  new_mesh->mpc->mpc_item[i] = mesh->mpc->mpc_item[i];
2807  new_mesh->mpc->mpc_dof[i] = mesh->mpc->mpc_dof[i];
2808  new_mesh->mpc->mpc_val[i] = mesh->mpc->mpc_val[i];
2809  }
2810  }
2811  }
2812  if (mesh->amp != NULL) {
2813  new_mesh->amp = (struct hecmwST_amplitude *)calloc(
2814  1, sizeof(struct hecmwST_amplitude));
2815  if (new_mesh->amp == NULL) HECMW_dlb_memory_exit("new_mesh: amp");
2816  new_mesh->amp->n_amp = mesh->amp->n_amp;
2817  if (new_mesh->amp->n_amp != 0) {
2819  (int *)calloc(mesh->amp->n_amp, sizeof(int));
2821  (int *)calloc(mesh->amp->n_amp, sizeof(int));
2823  (int *)calloc(mesh->amp->n_amp, sizeof(int));
2824  new_mesh->amp->amp_index =
2825  (int *)calloc(mesh->amp->n_amp + 1, sizeof(int));
2826  new_mesh->amp->amp_val = (double *)calloc(
2827  mesh->amp->amp_index[mesh->amp->n_amp], sizeof(double));
2828  new_mesh->amp->amp_table = (double *)calloc(
2829  mesh->amp->amp_index[mesh->amp->n_amp], sizeof(double));
2830  if ((new_mesh->amp->amp_type_definition == NULL) ||
2831  (new_mesh->amp->amp_type_time == NULL) ||
2832  (new_mesh->amp->amp_type_value == NULL) ||
2833  (new_mesh->amp->amp_index == NULL) ||
2834  (new_mesh->amp->amp_val == NULL) ||
2835  (new_mesh->amp->amp_table == NULL))
2836  HECMW_dlb_memory_exit("new_mesh: amp");
2837  for (i = 0; i < mesh->amp->n_amp; i++) {
2842  }
2843  for (i = 0; i < mesh->amp->n_amp + 1; i++)
2844  new_mesh->amp->amp_index[i] = mesh->amp->amp_index[i];
2845  for (i = 0; i < mesh->amp->amp_index[mesh->amp->n_amp]; i++) {
2846  new_mesh->amp->amp_val[i] = mesh->amp->amp_val[i];
2847  new_mesh->amp->amp_table[i] = mesh->amp->amp_table[i];
2848  }
2849  }
2850  }
2851  }
2853  (int *)calloc(new_mesh->nn_internal, sizeof(int));
2855  HECMW_dlb_memory_exit("node_internal_list");
2856  for (i = 0; i < new_mesh->nn_internal; i++)
2857  new_mesh->node_internal_list[i] = i;
2858  if (mynode == 0) {
2859  t2 = HECMW_Wtime();
2860  fprintf(stderr,
2861  "Finish migration now. The time cost for migration and generating "
2862  "new mesh is %lf\n",
2863  t2 - t1);
2864  }
2865  return;
2866 }
hecmwST_local_mesh::my_rank
int my_rank
Definition: hecmw_struct.h:212
hecmwST_local_mesh::mpc
struct hecmwST_mpc * mpc
Definition: hecmw_struct.h:247
hecmwST_result_data
Definition: hecmw_result.h:11
hecmwST_local_mesh::global_node_ID
int * global_node_ID
Definition: hecmw_struct.h:168
hecmwST_local_mesh::section_ID
int * section_ID
Definition: hecmw_struct.h:197
hecmwST_mpc
Definition: hecmw_struct.h:48
hecmwST_local_mesh::hecmw_n_file
int hecmw_n_file
Definition: hecmw_struct.h:155
hecmwST_local_mesh::shared_index
int * shared_index
Definition: hecmw_struct.h:221
if
if(!(yy_init))
Definition: hecmw_ablex.c:1823
hecmwST_local_mesh::errnof
int errnof
Definition: hecmw_struct.h:213
int2_whole_send_recv
void int2_whole_send_recv(int n1, int n2, int pesize, int *stack_import, int *stack_export, int *x, int *y, HECMW_Comm repart_comm, int my_rank)
Definition: hecmw_dlb_comm_util.c:519
hecmwST_amplitude::amp_val
double * amp_val
Definition: hecmw_struct.h:71
hecmwST_node_grp::n_grp
int n_grp
Definition: hecmw_struct.h:76
hecmwST_section::sect_R_index
int * sect_R_index
Definition: hecmw_struct.h:34
HECMW_dlb_memory_exit
void HECMW_dlb_memory_exit(char *var)
Definition: hecmw_dlb_mem_util.c:7
hecmwST_amplitude::amp_index
int * amp_index
Definition: hecmw_struct.h:70
HECMW_Allreduce
int HECMW_Allreduce(void *sendbuf, void *recvbuf, int count, HECMW_Datatype datatype, HECMW_Op op, HECMW_Comm comm)
Definition: hecmw_comm.c:364
hecmwST_material::n_mat_subitem
int n_mat_subitem
Definition: hecmw_struct.h:38
stack_whole_send_recv
int stack_whole_send_recv(int pesize, int *stack_import, int *stack_export, HECMW_Comm repart_comm, int my_rank)
Definition: hecmw_dlb_comm_util.c:173
HECMW_Send
int HECMW_Send(void *buffer, int count, HECMW_Datatype datatype, int dest, int tag, HECMW_Comm comm)
Definition: hecmw_comm.c:193
hecmwST_local_mesh::hecmw_flag_partdepth
int hecmw_flag_partdepth
Definition: hecmw_struct.h:146
hecmwST_material::n_mat_table
int n_mat_table
Definition: hecmw_struct.h:39
hecmwST_local_mesh::elem_node_item
int * elem_node_item
Definition: hecmw_struct.h:196
hecmwST_amplitude::amp_table
double * amp_table
Definition: hecmw_struct.h:72
hecmwST_local_mesh::gridfile
char gridfile[HECMW_FILENAME_LEN+1]
Definition: hecmw_struct.h:154
hecmwST_local_mesh::elem_mat_int_val
double * elem_mat_int_val
Definition: hecmw_struct.h:203
hecmwST_local_mesh::elem_group
struct hecmwST_elem_grp * elem_group
Definition: hecmw_struct.h:250
hecmwST_local_mesh::adapt_parent
int * adapt_parent
Definition: hecmw_struct.h:233
new_mesh
struct hecmwST_local_mesh * new_mesh
Definition: hecmw_repart.h:72
hecmwST_local_mesh::adapt_parent_type
int * adapt_parent_type
Definition: hecmw_struct.h:229
hecmwST_section::sect_I_index
int * sect_I_index
Definition: hecmw_struct.h:32
hecmwST_local_mesh::elem_type_item
int * elem_type_item
Definition: hecmw_struct.h:194
hecmwST_local_mesh::node_ID
int * node_ID
Definition: hecmw_struct.h:167
hecmwST_local_mesh::hecmw_flag_parttype
int hecmw_flag_parttype
Definition: hecmw_struct.h:142
mesh_migration_adapt
void mesh_migration_adapt(int mynode, int pesize, Result_part *result, int *vtxdist)
Definition: hecmw_dlb_elem_sr_adapt.c:50
hecmwST_mpc::mpc_item
int * mpc_item
Definition: hecmw_struct.h:51
hecmwST_local_mesh
Definition: hecmw_struct.h:139
_tmp_grp_inf
Definition: hecmw_repart.h:29
hecmwST_local_mesh::elem_ID
int * elem_ID
Definition: hecmw_struct.h:189
hecmwST_local_mesh::adapt_type
int * adapt_type
Definition: hecmw_struct.h:230
hecmwST_amplitude::amp_type_time
int * amp_type_time
Definition: hecmw_struct.h:64
hecmwST_amplitude::amp_type_value
int * amp_type_value
Definition: hecmw_struct.h:67
_tmp_grp_inf::item
int * item
Definition: hecmw_repart.h:31
hecmwST_local_mesh::elem_type
int * elem_type
Definition: hecmw_struct.h:191
hecmwST_mpc::mpc_val
double * mpc_val
Definition: hecmw_struct.h:53
hecmw_repart.h
hecmwST_local_mesh::n_elem
int n_elem
Definition: hecmw_struct.h:184
hecmwST_surf_grp::grp_index
int * grp_index
Definition: hecmw_struct.h:109
hecmwST_material::mat_subitem_index
int * mat_subitem_index
Definition: hecmw_struct.h:42
hecmwST_local_mesh::neighbor_pe
int * neighbor_pe
Definition: hecmw_struct.h:216
hecmwST_local_mesh::ne_internal
int ne_internal
Definition: hecmw_struct.h:186
hecmwST_local_mesh::export_item
int * export_item
Definition: hecmw_struct.h:220
int_whole_send_recv
void int_whole_send_recv(int n1, int n2, int pesize, int *stack_import, int *nod_import, int *stack_export, int *nod_export, int *x, int *y, HECMW_Comm repart_comm, int my_rank)
Definition: hecmw_dlb_comm_util.c:421
stack_part_send_recv
int stack_part_send_recv(int neibpetot, int *neibpe, int *stack_import, int *stack_export, HECMW_Comm repart_comm, int my_rank)
Definition: hecmw_dlb_comm_util.c:128
HECMW_Comm_dup
int HECMW_Comm_dup(HECMW_Comm comm, HECMW_Comm *new_comm)
Definition: hecmw_comm.c:56
hecmwST_local_mesh::node_dof_item
int * node_dof_item
Definition: hecmw_struct.h:175
whole_copy_array
void whole_copy_array(int *recv_elem_num, int *global_recv_elem_num, int mynode, int pesize, HECMW_Comm repart_comm)
Definition: hecmw_dlb_comm_util.c:100
hecmwST_mpc::mpc_dof
int * mpc_dof
Definition: hecmw_struct.h:52
hecmwST_local_mesh::export_index
int * export_index
Definition: hecmw_struct.h:219
hecmwST_material::mat_table_index
int * mat_table_index
Definition: hecmw_struct.h:43
int_part_send_recv
int int_part_send_recv(int n, int neibpetot, int *neibpe, int *stack_import, int *nod_import, int *stack_export, int *nod_export, int *x, HECMW_Comm repart_comm, int my_rank)
Definition: hecmw_dlb_comm_util.c:260
hecmwST_local_mesh::import_item
int * import_item
Definition: hecmw_struct.h:218
hecmwST_surf_grp::grp_item
int * grp_item
Definition: hecmw_struct.h:111
HECMW_SUM
#define HECMW_SUM
Definition: hecmw_config.h:58
data
struct hecmwST_result_data * data
Definition: neu_reporter.cpp:18
hecmwST_result_data::ne_component
int ne_component
Definition: hecmw_result.h:17
hecmwST_local_mesh::coarse_grid_level
int coarse_grid_level
Definition: hecmw_struct.h:225
hecmwST_local_mesh::n_node
int n_node
Definition: hecmw_struct.h:161
hecmwST_local_mesh::node
double * node
Definition: hecmw_struct.h:170
hecmwST_section::sect_opt
int * sect_opt
Definition: hecmw_struct.h:23
hecmwST_amplitude
Definition: hecmw_struct.h:57
hecmwST_elem_grp
Definition: hecmw_struct.h:92
mesh
struct hecmwST_local_mesh * mesh
Definition: hecmw_repart.h:71
hecmwST_local_mesh::PETOT
int PETOT
Definition: hecmw_struct.h:210
HECMW_Status
MPI_Status HECMW_Status
Definition: hecmw_config.h:36
hecmwST_local_mesh::when_i_was_refined_elem
int * when_i_was_refined_elem
Definition: hecmw_struct.h:228
hecmwST_local_mesh::node_dof_index
int * node_dof_index
Definition: hecmw_struct.h:174
hecmwST_local_mesh::n_elem_mat_ID
int n_elem_mat_ID
Definition: hecmw_struct.h:200
hecmwST_section::sect_mat_ID_index
int * sect_mat_ID_index
Definition: hecmw_struct.h:30
hecmwST_amplitude::n_amp
int n_amp
Definition: hecmw_struct.h:58
hecmwST_local_mesh::material
struct hecmwST_material * material
Definition: hecmw_struct.h:246
hecmwST_local_mesh::elem_mat_ID_index
int * elem_mat_ID_index
Definition: hecmw_struct.h:198
hecmwST_material::mat_temp
double * mat_temp
Definition: hecmw_struct.h:45
hecmwST_local_mesh::hecmw_flag_adapt
int hecmw_flag_adapt
Definition: hecmw_struct.h:140
hecmwST_mpc::n_mpc
int n_mpc
Definition: hecmw_struct.h:49
hecmwST_surf_grp::n_grp
int n_grp
Definition: hecmw_struct.h:107
double_whole_send_recv
void double_whole_send_recv(int n1, int n2, int pesize, int *stack_import, int *nod_import, int *stack_export, int *nod_export, double *x, double *y, HECMW_Comm repart_comm, int my_rank)
Definition: hecmw_dlb_comm_util.c:705
hecmwST_local_mesh::amp
struct hecmwST_amplitude * amp
Definition: hecmw_struct.h:248
hecmwST_local_mesh::zero_temp
double zero_temp
Definition: hecmw_struct.h:158
hecmwST_local_mesh::node_internal_list
int * node_internal_list
Definition: hecmw_struct.h:165
hecmwST_local_mesh::adapt_level
int * adapt_level
Definition: hecmw_struct.h:231
hecmwST_local_mesh::elem_node_index
int * elem_node_index
Definition: hecmw_struct.h:195
hecmwST_local_mesh::n_elem_type
int n_elem_type
Definition: hecmw_struct.h:192
hecmwST_section::sect_type
int * sect_type
Definition: hecmw_struct.h:17
HECMW_ANY_TAG
int HECMW_ANY_TAG
Definition: hecmw_dlb_comm_util.h:7
hecmwST_section
Definition: hecmw_struct.h:11
hecmwST_result_data::node_val_item
double * node_val_item
Definition: hecmw_result.h:25
hecmwST_section::n_sect
int n_sect
Definition: hecmw_struct.h:15
hecmwST_local_mesh::PEsmpTOT
int PEsmpTOT
Definition: hecmw_struct.h:211
hecmwST_elem_grp::grp_item
int * grp_item
Definition: hecmw_struct.h:96
hecmwST_node_grp::grp_index
int * grp_index
Definition: hecmw_struct.h:78
hecmwST_elem_grp::grp_name
char ** grp_name
Definition: hecmw_struct.h:94
hecmwST_section::sect_I_item
int * sect_I_item
Definition: hecmw_struct.h:33
hecmwST_local_mesh::hecmw_flag_initcon
int hecmw_flag_initcon
Definition: hecmw_struct.h:141
hecmwST_local_mesh::adapt_children_index
int * adapt_children_index
Definition: hecmw_struct.h:234
hecmwST_local_mesh::elem_type_index
int * elem_type_index
Definition: hecmw_struct.h:193
hecmwST_result_data::nn_dof
int * nn_dof
Definition: hecmw_result.h:19
hecmwST_result_data::nn_component
int nn_component
Definition: hecmw_result.h:16
hecmwST_result_data::elem_label
char ** elem_label
Definition: hecmw_result.h:23
hecmwST_node_grp::grp_item
int * grp_item
Definition: hecmw_struct.h:79
_result_partition_struct
Definition: hecmw_repart.h:52
repart_comm
int repart_comm
Definition: hecmw_repart.h:70
hecmwST_section::sect_mat_ID_item
int * sect_mat_ID_item
Definition: hecmw_struct.h:31
hecmwST_material::n_mat
int n_mat
Definition: hecmw_struct.h:36
hecmwST_node_grp
Definition: hecmw_struct.h:75
hecmwST_local_mesh::elem_mat_ID_item
int * elem_mat_ID_item
Definition: hecmw_struct.h:199
hecmwST_material::n_mat_item
int n_mat_item
Definition: hecmw_struct.h:37
new_data
struct hecmwST_result_data * new_data
Definition: hecmw_repart.h:76
hecmwST_local_mesh::when_i_was_refined_node
int * when_i_was_refined_node
Definition: hecmw_struct.h:227
hecmwST_material::mat_name
char ** mat_name
Definition: hecmw_struct.h:40
hecmwST_local_mesh::nn_internal
int nn_internal
Definition: hecmw_struct.h:164
hecmwST_local_mesh::n_adapt
int n_adapt
Definition: hecmw_struct.h:226
_tmp_grp_inf::num_of_item
int num_of_item
Definition: hecmw_repart.h:30
hecmwST_surf_grp
Definition: hecmw_struct.h:106
hecmwST_material::mat_item_index
int * mat_item_index
Definition: hecmw_struct.h:41
HECMW_Comm
MPI_Comm HECMW_Comm
Definition: hecmw_config.h:30
hecmwST_local_mesh::n_dof
int n_dof
Definition: hecmw_struct.h:171
double_part_send_recv
int double_part_send_recv(int n, int neibpetot, int *neibpe, int *stack_import, int *nod_import, int *stack_export, int *nod_export, double *x, HECMW_Comm repart_comm, int my_rank)
Definition: hecmw_dlb_comm_util.c:336
hecmwST_local_mesh::n_subdomain
int n_subdomain
Definition: hecmw_struct.h:214
hecmwST_result_data::node_label
char ** node_label
Definition: hecmw_result.h:22
hecmwST_local_mesh::import_index
int * import_index
Definition: hecmw_struct.h:217
hecmwST_local_mesh::hecmw_flag_version
int hecmw_flag_version
Definition: hecmw_struct.h:147
HECMW_INT
#define HECMW_INT
Definition: hecmw_config.h:48
HECMW_Wtime
double HECMW_Wtime(void)
Definition: hecmw_time.c:8
hecmwST_local_mesh::zero
int zero
Definition: hecmw_struct.h:208
hecmwST_local_mesh::files
char ** files
Definition: hecmw_struct.h:156
hecmwST_elem_grp::n_grp
int n_grp
Definition: hecmw_struct.h:93
HECMW_Recv
int HECMW_Recv(void *buffer, int count, HECMW_Datatype datatype, int source, int tag, HECMW_Comm comm, HECMW_Status *status)
Definition: hecmw_comm.c:235
hecmwST_material::mat_val
double * mat_val
Definition: hecmw_struct.h:44
hecmwST_local_mesh::section
struct hecmwST_section * section
Definition: hecmw_struct.h:245
hecmwST_elem_grp::grp_index
int * grp_index
Definition: hecmw_struct.h:95
hecmwST_node_grp::grp_name
char ** grp_name
Definition: hecmw_struct.h:77
int3_whole_send_recv
void int3_whole_send_recv(int n1, int n2, int pesize, int *stack_import, int *stack_export, int *x, int *y, HECMW_Comm repart_comm, int my_rank)
Definition: hecmw_dlb_comm_util.c:664
NULL
#define NULL
Definition: hecmw_io_nastran.c:30
hecmwST_local_mesh::adapt_children_item
int * adapt_children_item
Definition: hecmw_struct.h:235
hecmwST_material
Definition: hecmw_struct.h:35
hecmwST_section::sect_R_item
double * sect_R_item
Definition: hecmw_struct.h:35
hecmwST_local_mesh::header
char header[HECMW_HEADER_LEN+1]
Definition: hecmw_struct.h:157
hecmwST_mpc::mpc_index
int * mpc_index
Definition: hecmw_struct.h:50
hecmwST_surf_grp::grp_name
char ** grp_name
Definition: hecmw_struct.h:108
hecmwST_local_mesh::n_neighbor_pe
int n_neighbor_pe
Definition: hecmw_struct.h:215
HECMW_dlb_print_exit
void HECMW_dlb_print_exit(char *var)
Definition: hecmw_dlb_mem_util.c:16
hecmwST_local_mesh::global_elem_ID
int * global_elem_ID
Definition: hecmw_struct.h:190
hecmwST_local_mesh::surf_group
struct hecmwST_surf_grp * surf_group
Definition: hecmw_struct.h:251
hecmwST_local_mesh::n_dof_grp
int n_dof_grp
Definition: hecmw_struct.h:172
hecmwST_local_mesh::node_group
struct hecmwST_node_grp * node_group
Definition: hecmw_struct.h:249
_result_partition_struct::part
int * part
Definition: hecmw_repart.h:55
HECMW_Barrier
int HECMW_Barrier(HECMW_Comm comm)
Definition: hecmw_comm.c:95
double2_whole_send_recv
void double2_whole_send_recv(int n1, int n2, int pesize, int *stack_import, int *stack_export, double *x, double *y, HECMW_Comm repart_comm, int my_rank)
Definition: hecmw_dlb_comm_util.c:592
hecmwST_local_mesh::HECMW_COMM
HECMW_Comm HECMW_COMM
Definition: hecmw_struct.h:209
hecmwST_amplitude::amp_type_definition
int * amp_type_definition
Definition: hecmw_struct.h:61
hecmwST_result_data::ne_dof
int * ne_dof
Definition: hecmw_result.h:20
hecmwST_local_mesh::elem_internal_list
int * elem_internal_list
Definition: hecmw_struct.h:187
_result_partition_struct::t_node
int t_node
Definition: hecmw_repart.h:54