proteus  1.8.1
C/C++/Fortran libraries
DumpMesh.cpp
Go to the documentation of this file.
1 #include "DumpMesh.h"
2 #include <PCU.h>
3 
4 /* see flcbdfWrappersModule.cpp:
5  flcbdfWrappersConvertPUMIPartitionToPython
6  and cmeshToolsModule.cpp:
7  cmeshToolsBuildPythonMeshInterface */
8 
9 static void print_int_2d(int* a, int h, int w, FILE* f)
10 {
11  for (int i = 0; i < h; ++i) {
12  for (int j = 0; j < w; ++j)
13  fprintf(f, "%d ", a[i * w + j]);
14  fprintf(f, "\n");
15  }
16 }
17 
18 static void print_double_2d(double* a, int h, int w, FILE* f)
19 {
20  for (int i = 0; i < h; ++i) {
21  for (int j = 0; j < w; ++j)
22  fprintf(f, "%f ", a[i * w + j]);
23  fprintf(f, "\n");
24  }
25 }
26 
27 static void dump_proteus_mesh_header(Mesh* m, FILE* f)
28 {
29  fprintf(f, "nElements_global = %d\n", m->nElements_global);
30  fprintf(f, "nNodes_global = %d\n", m->nNodes_global);
31  fprintf(f, "nNodes_element = %d\n", m->nNodes_element);
32  fprintf(f, "nNodes_elementBoundary = %d\n", m->nNodes_elementBoundary);
33  fprintf(f, "nElementBoundaries_element = %d\n", m->nElementBoundaries_element);
34  fprintf(f, "nElementBoundaries_global = %d\n", m->nElementBoundaries_global);
35  fprintf(f, "nInteriorElementBoundaries_global = %d\n",
37  fprintf(f, "nExteriorElementBoundaries_global = %d\n",
39  fprintf(f, "max_nElements_node = %d\n", m->max_nElements_node);
40  fprintf(f, "nEdges_global = %d\n", m->nEdges_global);
41  fprintf(f, "max_nNodeNeighbors_node = %d\n", m->max_nNodeNeighbors_node);
42 }
43 
44 static void dump_proteus_subdomain(Mesh* m, FILE* f)
45 {
46  dump_proteus_mesh_header(m, f);
47  fprintf(f, "elementNodesArray:\n");
48  print_int_2d(m->elementNodesArray, m->nElements_global,
49  m->nNodes_element, f);
50  fprintf(f, "nodeElementOffsets:\n");
51  if (!m->nodeElementOffsets)
52  fprintf(f, " NULL\n");
53  else {
54  for (int i = 0; i <= m->nNodes_global; ++i)
55  fprintf(f, "%d\n", m->nodeElementOffsets[i]);
56  int total = m->nodeElementOffsets[m->nNodes_global];
57  assert(total == m->nElements_global * m->nNodes_element);
58  fprintf(f, "nodeElementsArray:\n");
59  for (int i = 0; i < total; ++i)
60  fprintf(f, "%d\n", m->nodeElementsArray[i]);
61  }
62  fprintf(f, "elementNeighborsArray:\n");
63  print_int_2d(m->elementNeighborsArray, m->nElements_global,
65  fprintf(f, "elementBoundariesArray:\n");
66  print_int_2d(m->elementBoundariesArray, m->nElements_global,
68  fprintf(f, "elementBoundaryNodesArray:\n");
71  fprintf(f, "elementBoundaryElementsArray:\n");
73  2, f);
74  fprintf(f, "elementBoundaryLocalElementBoundariesArray:\n");
77  fprintf(f, "interiorElementBoundariesArray:\n");
78  for (int i = 0; i < m->nInteriorElementBoundaries_global; ++i)
79  fprintf(f, "%d\n", m->interiorElementBoundariesArray[i]);
80  fprintf(f, "exteriorElementBoundariesArray:\n");
81  for (int i = 0; i < m->nExteriorElementBoundaries_global; ++i)
82  fprintf(f, "%d\n", m->exteriorElementBoundariesArray[i]);
83  fprintf(f, "edgeNodesArray:\n");
84  print_int_2d(m->edgeNodesArray, m->nEdges_global, 2, f);
85  fprintf(f, "nodeStarOffsets:\n");
86  if (!m->nodeStarOffsets)
87  fprintf(f, " NULL\n");
88  else {
89  for (int i = 0; i <= m->nNodes_global; ++i)
90  fprintf(f, "%d\n", m->nodeStarOffsets[i]);
91  int total = m->nodeStarOffsets[m->nNodes_global];
92  assert(total == m->nEdges_global * 2);
93  fprintf(f, "nodeStarArray:\n");
94  for (int i = 0; i < total; ++i)
95  fprintf(f, "%d\n", m->nodeStarArray[i]);
96  }
97  fprintf(f, "elementMaterialTypes:\n");
98  for (int i = 0; i < m->nElements_global; ++i)
99  fprintf(f, "%d\n", m->elementMaterialTypes[i]);
100  fprintf(f, "elementBoundaryMaterialTypes:\n");
101  for (int i = 0; i < m->nElementBoundaries_global; ++i)
102  fprintf(f, "%d\n", m->elementBoundaryMaterialTypes[i]);
103  fprintf(f, "nodeMaterialTypes:\n");
104  for (int i = 0; i < m->nNodes_global; ++i)
105  fprintf(f, "%d\n", m->nodeMaterialTypes[i]);
106  fprintf(f, "nodeArray:\n");
107  print_double_2d(m->nodeArray, m->nNodes_global, 3, f);
108  fprintf(f, "elementDiametersArray:\n");
109  for (int i = 0; i < m->nElements_global; ++i)
110  fprintf(f, "%f\n", m->elementDiametersArray[i]);
111  fprintf(f, "elementInnerDiametersArray:\n");
112  for (int i = 0; i < m->nElements_global; ++i)
113  fprintf(f, "%f\n", m->elementInnerDiametersArray[i]);
114  fprintf(f, "elementBoundaryDiametersArray:\n");
115  for (int i = 0; i < m->nElementBoundaries_global; ++i)
116  fprintf(f, "%f\n", m->elementBoundaryDiametersArray[i]);
117  fprintf(f, "elementBarycentersArray:\n");
118  print_double_2d(m->elementBarycentersArray, m->nElements_global, 3, f);
119  print_double_2d(m->elementBoundaryBarycentersArray,
121  fprintf(f, "nodeDiametersArray:\n");
122  for (int i = 0; i < m->nNodes_global; ++i)
123  fprintf(f, "%f\n", m->nodeDiametersArray[i]);
124  fprintf(f, "nodeSupportArray:\n");
125  for (int i = 0; i < m->nNodes_global; ++i)
126  fprintf(f, "%f\n", m->nodeSupportArray[i]);
127  fprintf(f, "h = %f\n", m->h);
128  fprintf(f, "hMin = %f\n", m->hMin);
129  fprintf(f, "sigmaMax = %f\n", m->sigmaMax);
130  fprintf(f, "volume = %f\n", m->volume);
131 }
132 
133 static void dump_proteus_parallel_arrays(Mesh* m, FILE* f)
134 {
135  int size = PCU_Comm_Peers();
136  fprintf(f, "elementOffsets_subdomain_owned:\n");
137  for (int i = 0; i <= size; ++i)
138  fprintf(f, "%d\n", m->elementOffsets_subdomain_owned[i]);
139  fprintf(f, "elementNumbering_subdomain2global:\n");
140  for (int i = 0; i < m->subdomainp->nElements_global; ++i)
141  fprintf(f, "%d\n", m->elementNumbering_subdomain2global[i]);
142  fprintf(f, "nodeOffsets_subdomain_owned:\n");
143  for (int i = 0; i <= size; ++i)
144  fprintf(f, "%d\n", m->nodeOffsets_subdomain_owned[i]);
145  fprintf(f, "nodeNumbering_subdomain2global:\n");
146  for (int i = 0; i < m->subdomainp->nNodes_global; ++i)
147  fprintf(f, "%d\n", m->nodeNumbering_subdomain2global[i]);
148  fprintf(f, "elementBoundaryOffsets_subdomain_owned:\n");
149  for (int i = 0; i <= size; ++i)
150  fprintf(f, "%d\n", m->elementBoundaryOffsets_subdomain_owned[i]);
151  fprintf(f, "elementBoundaryNumbering_subdomain2global:\n");
152  for (int i = 0; i < m->subdomainp->nElementBoundaries_global; ++i)
153  fprintf(f, "%d\n", m->elementBoundaryNumbering_subdomain2global[i]);
154  fprintf(f, "edgeOffsets_subdomain_owned:\n");
155  for (int i = 0; i <= size; ++i)
156  fprintf(f, "%d\n", m->edgeOffsets_subdomain_owned[i]);
157  fprintf(f, "edgeNumbering_subdomain2global:\n");
158  for (int i = 0; i < m->subdomainp->nEdges_global; ++i)
159  fprintf(f, "%d\n", m->edgeNumbering_subdomain2global[i]);
160 
161 }
162 
163 void dump_proteus_mesh(Mesh* m, FILE* f)
164 {
165  fprintf(stderr,"%d dumping mesh file\n", PCU_Comm_Self());
166  if (PCU_Comm_Peers() > 1) {
167  fprintf(f, "PARALLEL HEADER\n");
168  dump_proteus_mesh_header(m, f);
169  fprintf(f, "PARALLEL ARRAYS\n");
170  dump_proteus_parallel_arrays(m, f);
171  fprintf(f, "SUBDOMAIN\n");
172  dump_proteus_subdomain(m->subdomainp, f);
173  } else {
174  dump_proteus_subdomain(m, f);
175  }
176 }
Mesh::nodeDiametersArray
double * nodeDiametersArray
Definition: mesh.h:69
Mesh::elementBoundaryOffsets_subdomain_owned
int * elementBoundaryOffsets_subdomain_owned
Definition: mesh.h:80
Mesh::interiorElementBoundariesArray
int * interiorElementBoundariesArray
Definition: mesh.h:50
Mesh::elementBoundaryNodesArray
int * elementBoundaryNodesArray
Definition: mesh.h:47
Mesh::nodeArray
double * nodeArray
Definition: mesh.h:67
Mesh::nodeElementsArray
int * nodeElementsArray
Definition: mesh.h:43
Mesh::max_nElements_node
int max_nElements_node
Definition: mesh.h:38
w
#define w(x)
Definition: jf.h:22
Mesh::nodeElementOffsets
int * nodeElementOffsets
Definition: mesh.h:44
Mesh
Definition: mesh.h:28
Mesh::nEdges_global
int nEdges_global
Definition: mesh.h:39
Mesh::exteriorElementBoundariesArray
int * exteriorElementBoundariesArray
Definition: mesh.h:51
f
Double f
Definition: Headers.h:64
Mesh::nodeStarOffsets
int * nodeStarOffsets
Definition: mesh.h:54
Mesh::volume
double volume
Definition: mesh.h:70
Mesh::nNodes_elementBoundary
int nNodes_elementBoundary
Definition: mesh.h:33
Mesh::nElementBoundaries_global
int nElementBoundaries_global
Definition: mesh.h:35
Mesh::hMin
double hMin
Definition: mesh.h:70
Mesh::nInteriorElementBoundaries_global
int nInteriorElementBoundaries_global
Definition: mesh.h:36
Mesh::nodeOffsets_subdomain_owned
int * nodeOffsets_subdomain_owned
Definition: mesh.h:78
Mesh::nodeStarArray
int * nodeStarArray
Definition: mesh.h:53
Mesh::h
double h
Definition: mesh.h:70
Mesh::nNodes_element
int nNodes_element
Definition: mesh.h:32
Mesh::elementBoundaryLocalElementBoundariesArray
int * elementBoundaryLocalElementBoundariesArray
Definition: mesh.h:49
Mesh::elementBoundaryMaterialTypes
int * elementBoundaryMaterialTypes
Definition: mesh.h:56
Mesh::edgeNumbering_subdomain2global
int * edgeNumbering_subdomain2global
Definition: mesh.h:83
Mesh::elementInnerDiametersArray
double * elementInnerDiametersArray
Definition: mesh.h:67
Mesh::edgeOffsets_subdomain_owned
int * edgeOffsets_subdomain_owned
Definition: mesh.h:82
Mesh::elementNodesArray
int * elementNodesArray
Definition: mesh.h:42
Mesh::elementBoundaryBarycentersArray
double * elementBoundaryBarycentersArray
Definition: mesh.h:68
Mesh::elementBoundariesArray
int * elementBoundariesArray
Definition: mesh.h:46
Mesh::elementBarycentersArray
double * elementBarycentersArray
Definition: mesh.h:68
Mesh::nElements_global
int nElements_global
Definition: mesh.h:30
Mesh::nodeMaterialTypes
int * nodeMaterialTypes
Definition: mesh.h:57
Mesh::subdomainp
Mesh * subdomainp
Definition: mesh.h:84
Mesh::nElementBoundaries_element
int nElementBoundaries_element
Definition: mesh.h:34
Mesh::elementMaterialTypes
int * elementMaterialTypes
Definition: mesh.h:55
Mesh::elementOffsets_subdomain_owned
int * elementOffsets_subdomain_owned
Definition: mesh.h:76
Mesh::elementDiametersArray
double * elementDiametersArray
Definition: mesh.h:67
Mesh::elementBoundaryNumbering_subdomain2global
int * elementBoundaryNumbering_subdomain2global
Definition: mesh.h:81
Mesh::nNodes_global
int nNodes_global
Definition: mesh.h:31
Mesh::nExteriorElementBoundaries_global
int nExteriorElementBoundaries_global
Definition: mesh.h:37
DumpMesh.h
Mesh::max_nNodeNeighbors_node
int max_nNodeNeighbors_node
Definition: mesh.h:40
Mesh::nodeNumbering_subdomain2global
int * nodeNumbering_subdomain2global
Definition: mesh.h:79
Mesh::sigmaMax
double sigmaMax
Definition: mesh.h:70
Mesh::nodeSupportArray
double * nodeSupportArray
Definition: mesh.h:69
Mesh::elementNumbering_subdomain2global
int * elementNumbering_subdomain2global
Definition: mesh.h:77
dump_proteus_mesh
void dump_proteus_mesh(Mesh *m, FILE *f)
Definition: DumpMesh.cpp:163
Mesh::edgeNodesArray
int * edgeNodesArray
Definition: mesh.h:52
Mesh::elementBoundaryDiametersArray
double * elementBoundaryDiametersArray
Definition: mesh.h:67
Mesh::elementBoundaryElementsArray
int * elementBoundaryElementsArray
Definition: mesh.h:48
Mesh::elementNeighborsArray
int * elementNeighborsArray
Definition: mesh.h:45