gmsh-TingyuanDoc  0.1
An Open-Source Timing-driven Analytical Mixed-size FPGA Placer
CGNSZone.cpp
Go to the documentation of this file.
1 // Gmsh - Copyright (C) 1997-2022 C. Geuzaine, J.-F. Remacle
2 //
3 // See the LICENSE.txt file in the Gmsh root directory for license information.
4 // Please report all issues on https://gitlab.onelab.info/gmsh/gmsh/issues.
5 //
6 // Contributor(s):
7 // Thomas Toulorge
8 
9 #include "Context.h"
10 #include "GmshMessage.h"
11 #include "MVertex.h"
12 #include "affineTransformation.h"
13 #include "CGNSCommon.h"
14 #include "CGNSConventions.h"
15 #include "CGNSRead.h"
16 #include "CGNSZone.h"
17 
18 #if defined(HAVE_LIBCGNS)
19 
20 CGNSZone::CGNSZone(int fileIndex, int baseIndex, int zoneIndex,
21  CGNS_ENUMT(ZoneType_t) type, int meshDim, cgsize_t startNode,
22  const Family2EltNodeTransfo &allEltNodeTransfo, int &err)
23  : fileIndex_(fileIndex), baseIndex_(baseIndex), meshDim_(meshDim),
24  zoneIndex_(zoneIndex), type_(type), startNode_(startNode),
25  eltNodeTransfo_(nullptr), nbPerConnect_(0)
26 {
27  int cgnsErr;
28 
29  // read zone name & size
30  char name[CGNS_MAX_STR_LEN];
31  cgnsErr = cg_zone_read(fileIndex, baseIndex, zoneIndex, name, size_);
32  if(cgnsErr != CG_OK) err = cgnsError(__FILE__, __LINE__, fileIndex);
33  name_ = std::string(name);
34 
35  // read family name and retrieve element node tranformations (CPEX0045)
36  cgnsErr = cg_goto(fileIndex, baseIndex, "Zone_t", zoneIndex, "end");
37  if(cgnsErr != CG_OK) err = cgnsError(__FILE__, __LINE__, fileIndex);
38  char famName[CGNS_MAX_STR_LEN];
39  cgnsErr = cg_famname_read(famName);
40  if(cgnsErr != CG_NODE_NOT_FOUND) {
41  if(cgnsErr == CG_OK) {
42  auto it = allEltNodeTransfo.find(std::string(famName));
43  if(it != allEltNodeTransfo.end()) eltNodeTransfo_ = &(it->second);
44  }
45  else
46  err = cgnsError(__FILE__, __LINE__, fileIndex);
47  }
48 
49  err = 1;
50 }
51 
52 // read a boundary condition in a zone
53 int CGNSZone::readBoundaryCondition(int iZoneBC,
54  const std::vector<CGNSZone *> &allZones,
55  std::vector<std::string> &allGeomName)
56 {
57  int cgnsErr;
58 
59  // read general information on boundary condition
60  char rawBCName[CGNS_MAX_STR_LEN];
61  CGNS_ENUMT(BCType_t) bcType;
62  CGNS_ENUMT(PointSetType_t) ptSetType;
63  cgsize_t nbVal, normalSize;
64  CGNS_ENUMT(DataType_t) normalType;
65  int nbDataSet;
66  int normalIndex;
67  cgnsErr = cg_boco_info(fileIndex(), baseIndex(), index(), iZoneBC, rawBCName,
68  &bcType, &ptSetType, &nbVal, &normalIndex, &normalSize,
69  &normalType, &nbDataSet);
70  if(cgnsErr != CG_OK) return cgnsError(__FILE__, __LINE__, fileIndex());
71 
72  // get BC name from family linked to BC, use original BC name if not present,
73  // then retrieve BC index
74  std::string geomName;
75  cgnsErr = cg_goto(fileIndex(), baseIndex(), "Zone_t", index(), "ZoneBC_t", 1,
76  "BC_t", iZoneBC, "end");
77  if(cgnsErr != CG_OK) return cgnsError(__FILE__, __LINE__, fileIndex());
78  char rawFamilyName[CGNS_MAX_STR_LEN];
79  cgnsErr = cg_famname_read(rawFamilyName);
80  if(cgnsErr != CG_NODE_NOT_FOUND) {
81  if(cgnsErr == CG_OK)
82  geomName = std::string(rawFamilyName);
83  else
84  return cgnsError(__FILE__, __LINE__, fileIndex());
85  }
86  else
87  geomName = std::string(rawBCName);
88  const int indGeom = nameIndex(geomName, allGeomName);
89 
90  // read location of bnd. condition (type of mesh entity on which it applies)
91  CGNS_ENUMT(GridLocation_t) location;
92  cgnsErr = cg_boco_gridlocation_read(fileIndex(), baseIndex(), index(),
93  iZoneBC, &location);
94  if(cgnsErr != CG_OK) return cgnsError(__FILE__, __LINE__, fileIndex());
95 
96  // check that boundary condition is imposed on face elements
97  if((type() == CGNS_ENUMV(Unstructured)) && (meshDim() == 2) &&
98  (location != CGNS_ENUMV(CellCenter)) &&
99  (location != CGNS_ENUMV(EdgeCenter))) {
100  Msg::Warning("Boundary condition %s is specified on %s instead of "
101  "CellCenter/EdgeCenter in a 2D zone, skipping",
102  geomName.c_str(), cg_GridLocationName(location));
103  return 1;
104  }
105  else if((type() == CGNS_ENUMV(Unstructured)) && (meshDim() == 3) &&
106  (location != CGNS_ENUMV(CellCenter)) &&
107  (location != CGNS_ENUMV(FaceCenter))) {
108  Msg::Warning("Boundary condition %s is specified on %s instead of "
109  "CellCenter/FaceCenter in a 3D zone, skipping",
110  geomName.c_str(), cg_GridLocationName(location));
111  return 1;
112  }
113 
114  // read and store elements on which the BC is imposed
115  std::vector<cgsize_t> bcElt;
116  switch(ptSetType) {
117  case CGNS_ENUMV(ElementRange):
118  case CGNS_ENUMV(PointRange):
119  readBoundaryConditionRange(iZoneBC, bcElt);
120  break;
121  case CGNS_ENUMV(ElementList):
122  case CGNS_ENUMV(PointList):
123  readBoundaryConditionList(iZoneBC, nbVal, bcElt);
124  break;
125  default:
126  Msg::Error("Wrong point set type %s is for boundary condition %s",
127  cg_PointSetTypeName(ptSetType), geomName.c_str());
128  return 0;
129  break;
130  }
131  for(std::size_t iElt = 0; iElt < bcElt.size(); iElt++) {
132  elt2Geom()[bcElt[iElt]] = indGeom;
133  }
134 
135  return 1;
136 }
137 
138 int CGNSZone::readVertices(int dim, double scale,
139  std::vector<CGNSZone *> &allZones,
140  std::vector<MVertex *> &zoneVert)
141 {
142  int cgnsErr;
143 
144  // read dimension of coordinates and check consistency with base node
145  int dimZone;
146  cgnsErr = cg_ncoords(fileIndex(), baseIndex(), index(), &dimZone);
147  if(cgnsErr != CG_OK) return cgnsError(__FILE__, __LINE__, fileIndex());
148  if(dimZone > dim) {
149  Msg::Warning("%i coordinates in CGNS zone %i, while base has dimension %i,"
150  " discarding upper dimensions",
151  dimZone, index(), dim);
152  }
153  else if(dimZone < dim) {
154  Msg::Error("%i coordinates in CGNS zone %i, while base has dimension %i",
155  dimZone, index(), dim);
156  return 0;
157  }
158 
159  // read vertex coordinates
160  std::vector<double> xyz[3];
161  for(int iXYZ = 0; iXYZ < dim; iXYZ++) {
162  char xyzName[CGNS_MAX_STR_LEN];
163  CGNS_ENUMT(DataType_t) dataType;
164  cgnsErr = cg_coord_info(fileIndex(), baseIndex(), index(), iXYZ + 1,
165  &dataType, xyzName);
166  if(cgnsErr != CG_OK) return cgnsError(__FILE__, __LINE__, fileIndex());
167  const cgsize_t startInd[3] = {1, 1, 1};
168  xyz[iXYZ].resize(nbNode());
169  cgnsErr =
170  cg_coord_read(fileIndex(), baseIndex(), index(), xyzName,
171  CGNS_ENUMV(RealDouble), startInd, size(), xyz[iXYZ].data());
172  if(cgnsErr != CG_OK) return cgnsError(__FILE__, __LINE__, fileIndex());
173  }
174 
175  // create vertices
176  zoneVert.reserve(nbNode());
177  for(int i = 0; i < nbNode(); i++) {
178  const double x = xyz[0][i] * scale;
179  const double y = (dim > 1) ? xyz[1][i] * scale : 0.;
180  const double z = (dim > 2) ? xyz[2][i] * scale : 0.;
181  zoneVert.push_back(new MVertex(x, y, z));
182  }
183 
184  return 1;
185 }
186 
187 int CGNSZone::readConnectivities(const std::map<std::string, int> &name2Zone,
188  std::vector<CGNSZone *> &allZones)
189 {
190  int cgnsErr;
191 
192  // read number of connectivities
193  int nbConnect;
194  cgnsErr = cg_nconns(fileIndex(), baseIndex(), index(), &nbConnect);
195  if(cgnsErr != CG_OK) return cgnsError(__FILE__, __LINE__, fileIndex());
196 
197  // reserve memory for data containers
198  perTransfo_.reserve(nbConnect);
199  slaveNode_.reserve(nbConnect);
200  slaveVert_.reserve(nbConnect);
201  masterZone_.reserve(nbConnect);
202  masterNode_.reserve(nbConnect);
203  masterVert_.reserve(nbConnect);
204 
205  for(int iConnect = 1; iConnect <= nbConnect; iConnect++) {
206  // read connection info
207  char connectName[CGNS_MAX_STR_LEN], donorName[CGNS_MAX_STR_LEN];
208  CGNS_ENUMT(GridLocation_t) location;
209  CGNS_ENUMT(GridConnectivityType_t) connectType;
210  CGNS_ENUMT(PointSetType_t) ptSetType, ptSetTypeDonor;
211  cgsize_t connectSize, connectSizeDonor;
212  CGNS_ENUMT(ZoneType_t) zoneTypeDonor;
213  CGNS_ENUMT(DataType_t) dataTypeDonor;
214  cgnsErr = cg_conn_info(fileIndex(), baseIndex(), index(), iConnect,
215  connectName, &location, &connectType, &ptSetType,
216  &connectSize, donorName, &zoneTypeDonor,
217  &ptSetTypeDonor, &dataTypeDonor, &connectSizeDonor);
218  if(cgnsErr != CG_OK) return cgnsError(__FILE__, __LINE__, fileIndex());
219 
220  // read periodic info, skip if not periodic
221  float rotCenter[3], rotAngle[3], trans[3];
222  cgnsErr = cg_conn_periodic_read(fileIndex(), baseIndex(), index(), iConnect,
223  rotCenter, rotAngle, trans);
224  if(cgnsErr == CG_NODE_NOT_FOUND) continue;
225  if(cgnsErr != CG_OK) return cgnsError(__FILE__, __LINE__, fileIndex());
226 
227  // invert transformation because CGNS transformation is from current zone
228  // to donor zone, while Gmsh transformation is from master to slave
229  for(int d = 0; d < 3; d++) {
230  rotAngle[d] = -rotAngle[d];
231  trans[d] = -trans[d];
232  }
233 
234  // check if connection type is OK
235  if(connectType != CGNS_ENUMV(Abutting1to1)) {
236  Msg::Warning("Non-conformal connection not supported in CGNS reader");
237  continue;
238  }
239  if(location != CGNS_ENUMV(Vertex)) {
240  Msg::Warning("Only vertex connections are supported in CGNS reader");
241  continue;
242  }
243 
244  // get and check data on master zone
245  const std::string masterName(donorName);
246  const auto itMasterName = name2Zone.find(masterName);
247  if(itMasterName == name2Zone.end()) {
248  Msg::Error("Zone name '%s' in not found in connection %i of zone %i",
249  masterName.c_str(), iConnect, index());
250  return 0;
251  }
252  const int masterZoneIndex = itMasterName->second;
253  CGNSZone *mZone = allZones[masterZoneIndex];
254  if(mZone->type() != zoneTypeDonor) {
255  Msg::Error("Inconsistent type for zone '%s' in connection %i of zone %i",
256  masterName.c_str(), iConnect, index());
257  return 0;
258  }
259 
260  // read connectivity data
261  std::vector<cgsize_t> slaveData(indexDataSize(connectSize));
262  std::vector<cgsize_t> masterData(mZone->indexDataSize(connectSizeDonor));
263  cgnsErr = cg_conn_read(fileIndex(), baseIndex(), index(), iConnect,
264  slaveData.data(), dataTypeDonor, masterData.data());
265  if(cgnsErr != CG_OK) return cgnsError(__FILE__, __LINE__, fileIndex());
266 
267  // get slave and master nodes
268  std::vector<cgsize_t> sNode, mNode;
269  if(ptSetType == CGNS_ENUMV(PointRange))
270  nodeFromRange(slaveData, sNode);
271  else if(ptSetType == CGNS_ENUMV(PointList))
272  nodeFromList(slaveData, sNode);
273  if(ptSetTypeDonor != CGNS_ENUMV(PointListDonor)) {
274  Msg::Error("Only PointListDonor sets are supported for donnor points for "
275  "general connections in CGNS reader");
276  return 0;
277  }
278  mZone->nodeFromList(masterData, mNode);
279 
280  // add periodic connection
281  nbPerConnect_++;
282  perTransfo_.push_back(std::vector<double>());
283  computeAffineTransformation(rotCenter, rotAngle, trans, perTransfo_.back());
284  masterZone_.push_back(masterZoneIndex);
285  slaveNode_.push_back(sNode);
286  slaveVert_.push_back(std::vector<MVertex *>());
287  masterNode_.push_back(mNode);
288  masterVert_.push_back(std::vector<MVertex *>());
289  }
290 
291  return 1;
292 }
293 
294 int CGNSZone::readMesh(int dim, double scale, std::vector<CGNSZone *> &allZones,
295  std::vector<MVertex *> &allVert,
296  std::map<int, std::vector<MElement *> > *allElt,
297  std::vector<MVertex *> &zoneVert,
298  std::vector<MElement *> &zoneElt,
299  std::vector<std::string> &allGeomName)
300 {
301  // read boundary conditions for classification of mesh on geometry
302  if(CTX::instance()->mesh.cgnsImportIgnoreBC == 0) {
303  int cgnsErr;
304  int nbZoneBC;
305  cgnsErr = cg_nbocos(fileIndex(), baseIndex(), index(), &nbZoneBC);
306  if(cgnsErr != CG_OK) return cgnsError(__FILE__, __LINE__, fileIndex());
307  for(int iZoneBC = 1; iZoneBC <= nbZoneBC; iZoneBC++) {
308  int errBC = readBoundaryCondition(iZoneBC, allZones, allGeomName);
309  if(errBC == 0) return 0;
310  }
311  }
312 
313  // read and create vertices
314  int errVert = readVertices(dim, scale, allZones, zoneVert);
315  if(errVert == 0) return 0;
316  allVert.insert(allVert.end(), zoneVert.begin(), zoneVert.end());
317 
318  // read and create elements
319  int err = readElements(allVert, allElt, zoneElt, allGeomName);
320  if(err == 0) return 0;
321 
322  // cleanup unncessary memory
323  elt2Geom().clear();
324 
325  return 1;
326 }
327 
328 void CGNSZone::setPeriodicVertices(const std::vector<CGNSZone *> &allZones,
329  const std::vector<MVertex *> &allVert)
330 {
331  for(int iPer = 0; iPer < nbPerConnect(); iPer++) {
332  const std::vector<cgsize_t> &sNode = slaveNode(iPer);
333  const std::vector<cgsize_t> &mNode = masterNode(iPer);
334  std::vector<MVertex *> &sVert = slaveVert(iPer);
335  std::vector<MVertex *> &mVert = masterVert(iPer);
336  CGNSZone *mZone = allZones[masterZone(iPer)];
337  for(std::size_t iN = 0; iN < sNode.size(); iN++) {
338  const cgsize_t sInd = startNode() + sNode[iN];
339  const cgsize_t mInd = mZone->startNode() + mNode[iN];
340  sVert.push_back(allVert[sInd]);
341  mVert.push_back(allVert[mInd]);
342  }
343  }
344 }
345 
346 int CGNSZone::readBoundaryConditionRange(int iZoneBC,
347  std::vector<cgsize_t> &bcElt)
348 {
349  int cgnsErr;
350 
351  std::vector<cgsize_t> bcData(indexDataSize(2));
352  cgnsErr = cg_boco_read(fileIndex(), baseIndex(), index(), iZoneBC,
353  bcData.data(), nullptr);
354  if(cgnsErr != CG_OK) return cgnsError(__FILE__, __LINE__, fileIndex());
355 
356  // get list of elements from range data
357  eltFromRange(bcData, bcElt);
358 
359  return 1;
360 }
361 
362 int CGNSZone::readBoundaryConditionList(int iZoneBC, cgsize_t nbVal,
363  std::vector<cgsize_t> &bcElt)
364 {
365  int cgnsErr;
366 
367  // read data
368  std::vector<cgsize_t> bcData(indexDataSize(nbVal));
369  cgnsErr = cg_boco_read(fileIndex(), baseIndex(), index(), iZoneBC,
370  bcData.data(), nullptr);
371  if(cgnsErr != CG_OK) return cgnsError(__FILE__, __LINE__, fileIndex());
372 
373  // get list of elements from list data
374  eltFromList(bcData, bcElt);
375 
376  return 1;
377 }
378 
379 #endif // HAVE_LIBCGNS
CGNSRead.h
MVertex
Definition: MVertex.h:24
Msg::Warning
static void Warning(const char *fmt,...)
Definition: GmshMessage.cpp:543
Msg::Error
static void Error(const char *fmt,...)
Definition: GmshMessage.cpp:482
Vertex
Definition: Geo.h:29
CGNSZone.h
GmshMessage.h
CTX::instance
static CTX * instance()
Definition: Context.cpp:122
affineTransformation.h
MVertex.h
computeAffineTransformation
bool computeAffineTransformation(const FLOAT *rc, const FLOAT *ra, const FLOAT *tr, std::vector< double > &tfo)
Definition: affineTransformation.cpp:14
CGNSConventions.h
CGNSCommon.h
Context.h
z
const double z
Definition: GaussQuadratureQuad.cpp:56
GModel::mesh
int mesh(int dimension)
Definition: GModel.cpp:1066
scale
static void scale(std::vector< double > &x, double s)
Definition: ConjugateGradients.cpp:21