FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
output_ucd_c.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 
10  int i, j, k;
11  FILE *fp;
12  int id_elem;
13  int flag_hit;
14  int tmp_count, n_elem_type;
15  int tn_component, ne_internal, tmp_int;
16 
17  if ((fp = fopen(outfile, "w")) == NULL) {
18  fprintf(stderr, "output file error: %s\n", outfile);
19  exit(1);
20  }
21 
22  tn_component = 0;
23  for (i = 0; i < new_data->nn_component; i++)
24  tn_component += new_data->nn_dof[i];
25  ne_internal = 0;
26  for (i = 0; i < new_mesh->ne_internal; i++) {
28  ne_internal++;
29  }
30 
31  fprintf(fp, "%d %d %d %d %d\n", new_mesh->n_node, ne_internal, tn_component,
32  0, 0);
33  for (i = 0; i < new_mesh->n_node; i++)
34  fprintf(fp, "%d %lf %lf %lf\n", i + 1, new_mesh->node[i * 3],
35  new_mesh->node[i * 3 + 1], new_mesh->node[i * 3 + 2]);
36  tmp_count = 0;
37  for (i = 0; i < new_mesh->ne_internal; i++) {
38  if (new_mesh->adapt_type[new_mesh->elem_internal_list[i] - 1] == 0) {
39  tmp_count++;
40  fprintf(fp, "%d 1 ", tmp_count);
41  tmp_int = new_mesh->elem_internal_list[i] - 1;
42  if ((new_mesh->elem_node_index[tmp_int + 1] -
43  new_mesh->elem_node_index[tmp_int]) == 4)
44  fprintf(fp, " tet ");
45  else
46  fprintf(fp, " prism ");
47  for (j = new_mesh->elem_node_index[tmp_int];
48  j < new_mesh->elem_node_index[tmp_int + 1]; j++)
49  fprintf(fp, "%d ", new_mesh->elem_node_item[j]);
50  fprintf(fp, "\n");
51  }
52  }
53  fprintf(fp, "%d ", new_data->nn_component);
54  for (i = 0; i < new_data->nn_component; i++)
55  fprintf(fp, "%d ", new_data->nn_dof[i]);
56  fprintf(fp, "\n");
57  for (i = 0; i < new_data->nn_component; i++)
58  fprintf(fp, "%s,\n", new_data->node_label[i]);
59  for (i = 0; i < new_mesh->n_node; i++) {
60  fprintf(fp, "%d ", i + 1);
61  for (j = 0; j < tn_component; j++)
62  fprintf(fp, "%lf ", new_data->node_val_item[i * tn_component + j]);
63  fprintf(fp, "\n");
64  }
65  fclose(fp);
66  return;
67 }
68 
69 /*
70 void write_one_mesh_display(char *outfile, struct hecmwST_local_mesh *new_mesh,
71  struct hecmwST_result_data *new_data,
72  MPI_Comm VIS_COMM, int mynode, int pesize) {
73  int i, j, k;
74  FILE *fp;
75  int id_elem;
76  int flag_hit;
77  int tmp_count, n_elem_type;
78  int tn_component, ne_internal, tmp_int;
79  int n_node, n_elem;
80  int nvtxs, *new_vtxdist, tmp_sum, tmp_nvtxs, *tmp_data;
81 
82  if (mynode == 0) {
83  if ((fp = fopen(outfile, "w")) == NULL) {
84  fprintf(stderr, "output file error: %s\n", outfile);
85  exit(1);
86  }
87  new_vtxdist = (int *)calloc(pesize + 1, sizeof(int));
88  nvtxs = new_mesh->nn_internal;
89  if (mynode == 0) {
90  new_vtxdist[0] = 0;
91  new_vtxdist[1] = nvtxs;
92  tmp_sum = nvtxs;
93  for (i = 1; i < pesize; i++) {
94  MPI_Recv(&tmp_nvtxs, 1, MPI_INT, i, MPI_ANY_TAG, VIS_COMM, &stat);
95  tmp_sum += tmp_nvtxs;
96  new_vtxdist[i + 1] = tmp_sum;
97  }
98  for (i = 1; i < pesize; i++)
99  MPI_Send(new_vtxdist, pesize + 1, MPI_INT, i, 0, VIS_COMM);
100  } else {
101  MPI_Send(&nvtxs, 1, MPI_INT, 0, 0, VIS_COMM);
102  MPI_Recv(new_vtxdist, pesize + 1, MPI_INT, 0, MPI_ANY_TAG, VIS_COMM,
103  &stat);
104  }
105 
106  tn_component = 0;
107  for (i = 0; i < new_data->nn_component; i++)
108  tn_component += new_data->nn_dof[i];
109  }
110  ne_internal = 0;
111  for (i = 0; i < new_mesh->ne_internal; i++) {
112  if (new_mesh->adapt_type[new_mesh->elem_internal_list[i] - 1] == 0)
113  ne_internal++;
114  }
115  if (pesize > 1)
116  MPI_Allreduce(&new_mesh->nn_internal, &n_node, 1, MPI_INT, MPI_SUM,
117  VIS_COMM);
118  else
119  n_node = mesh->nn_internal;
120  if (pesize > 1)
121  MPI_Allreduce(&ne_internal, &n_elem, 1, MPI_INT, MPI_SUM, VIS_COMM);
122  else
123  n_elem = ne_internal;
124 
125  if (mynode == 0) {
126  fprintf(fp, "%d %d %d %d %d\n", n_node, n_elem, tn_component, 0, 0);
127  for (i = 0; i < new_mesh->nn_internal; i++)
128  fprintf(fp, "%d %lf %lf %lf\n", i + 1, new_mesh->node[i * 3],
129  new_mesh->node[i * 3 + 1], new_mesh->node[i * 3 + 2]);
130  for (i = 1; i < pesize; i++) {
131  tmp_int = new_vtxdist[i] - new_vtxdist[i - 1];
132  tmp_data = (double *)calloc(tmp_int * 3, sizeof(double));
133  if (tmp_data == NULL) memory_exit("tmp_data");
134  MPI_Recv(tmp_data, tmp_int * 3, MPI_DOUBLE, i, MPI_ANY_TAG, VIS_COMM,
135  &stat);
136  for (j = 0; j < tmp_int; j++)
137  fprintf(fp, "%d %lf %lf %lf\n", new_vtxdist[i - 1] + j + 1,
138  tmp_data[j * 3], tmp_data[j * 3 + 1], tmp_data[j * 3 + 2]);
139  free(tmp_data);
140  }
141  } else if (mynode != 0) {
142  tmp_data = (double *)calloc(new_mesh->nn_internal * 3, sizeof(double));
143  if (tmp_data == NULL) memory_exit("tmp_data");
144  for (j = 0; j < new_mesh->nn_internal; j++) {
145  for (k = 0; k < 3; k++)
146  tmp_data[j*3+k]=new_mesh->node[j*3+k]);
147  }
148  MPI_Send(tmp_data, new_mesh->nn_internal * 3, MPI_DOUBLE, 0, 0, VIS_COMM);
149  free(tmp_data);
150  }
151 
152  tmp_count = 0;
153  for (i = 0; i < new_mesh->ne_internal; i++) {
154  if (new_mesh->adapt_type[new_mesh->elem_internal_list[i] - 1] == 0) {
155  tmp_count++;
156  fprintf(fp, "%d 1 ", tmp_count);
157  tmp_int = new_mesh->elem_internal_list[i] - 1;
158  if ((new_mesh->elem_node_index[tmp_int + 1] -
159  new_mesh->elem_node_index[tmp_int]) == 4)
160  fprintf(fp, " tet ");
161  else
162  fprintf(fp, " prism ");
163  for (j = new_mesh->elem_node_index[tmp_int];
164  j < new_mesh->elem_node_index[tmp_int + 1]; j++)
165  fprintf(fp, "%d ", new_mesh->elem_node_item[j]);
166  fprintf(fp, "\n");
167  }
168  }
169  fprintf(fp, "%d ", new_data->nn_component);
170  for (i = 0; i < new_data->nn_component; i++)
171  fprintf(fp, "%d ", new_data->nn_dof[i]);
172  fprintf(fp, "\n");
173  for (i = 0; i < new_data->nn_component; i++)
174  fprintf(fp, "%s,\n", new_data->node_label[i]);
175  for (i = 0; i < new_mesh->n_node; i++) {
176  fprintf(fp, "%d ", i + 1);
177  for (j = 0; j < tn_component; j++)
178  fprintf(fp, "%lf ", new_data->node_val_item[i * tn_component + j]);
179  fprintf(fp, "\n");
180  }
181  fclose(fp);
182  return;
183 }
184 */
hecmwST_result_data
Definition: hecmw_result.h:11
new_mesh
struct hecmwST_local_mesh * new_mesh
Definition: hecmw_repart.h:72
hecmwST_local_mesh::elem_node_item
int * elem_node_item
Definition: hecmw_struct.h:196
hecmwST_local_mesh
Definition: hecmw_struct.h:139
hecmwST_local_mesh::adapt_type
int * adapt_type
Definition: hecmw_struct.h:230
hecmw_repart.h
hecmwST_local_mesh::ne_internal
int ne_internal
Definition: hecmw_struct.h:186
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_local_mesh::elem_node_index
int * elem_node_index
Definition: hecmw_struct.h:195
new_data
struct hecmwST_result_data * new_data
Definition: hecmw_repart.h:76
hecmwST_result_data::node_val_item
double * node_val_item
Definition: hecmw_result.h:25
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
write_dist_mesh_display
void write_dist_mesh_display(char *outfile, struct hecmwST_local_mesh *new_mesh, struct hecmwST_result_data *new_data)
Definition: output_ucd_c.c:8
hecmwST_result_data::node_label
char ** node_label
Definition: hecmw_result.h:22
NULL
#define NULL
Definition: hecmw_io_nastran.c:30
hecmwST_local_mesh::elem_internal_list
int * elem_internal_list
Definition: hecmw_struct.h:187