1 /*******************************************************************/
2 /*                               XDMF                              */
3 /*                   eXtensible Data Model and Format              */
4 /*                                                                 */
5 /*  Id : Id  */
6 /*  Date : $Date$ */
7 /*  Version : $Revision$ */
8 /*                                                                 */
9 /*  Author:                                                        */
10 /*     Kenneth Leiter                                              */
11 /*     kenneth.leiter@arl.army.mil                                 */
12 /*     US Army Research Laboratory                                 */
13 /*     Aberdeen Proving Ground, MD                                 */
14 /*                                                                 */
15 /*     Copyright @ 2009 US Army Research Laboratory                */
16 /*     All Rights Reserved                                         */
17 /*     See Copyright.txt or http://www.arl.hpc.mil/ice for details */
18 /*                                                                 */
19 /*     This software is distributed WITHOUT ANY WARRANTY; without  */
20 /*     even the implied warranty of MERCHANTABILITY or FITNESS     */
21 /*     FOR A PARTICULAR PURPOSE.  See the above copyright notice   */
22 /*     for more information.                                       */
23 /*                                                                 */
24 /*******************************************************************/
25 
26 #include <Xdmf.h>
27 #include <XdmfSet.h>
28 
29 #include <stdio.h>
30 #include <sstream>
31 #include <map>
32 #include <stack>
33 #include <vector>
34 
35 #include "XdmfFortran.h"
36 
37 // This works with g77. Different compilers require different
38 // name mangling.
39 #if !defined(WIN32)
40 #define XdmfInit xdmfinit_
41 #define XdmfSetTime xdmfsettime_
42 #define XdmfAddCollection xdmfaddcollection_
43 #define XdmfCloseCollection xdmfclosecollection_
44 #define XdmfSetGridTopology xdmfsetgridtopology_
45 #define XdmfSetGridTopologyFromShape xdmfsetgridtopologyfromshape_
46 #define XdmfSetGridGeometry xdmfsetgridgeometry_
47 #define XdmfAddGridAttribute xdmfaddgridattribute_
48 #define XdmfAddGridAttributeFromShape xdmfaddgridattributefromshape_
49 #define XdmfAddGridInformation xdmfaddgridinformation_
50 #define XdmfAddCollectionAttribute xdmfaddcollectionattribute_
51 #define XdmfAddCollectionInformation xdmfaddcollectioninformation_
52 #define XdmfAddArray xdmfaddarray_
53 #define XdmfReadFile xdmfreadfile_
54 #define XdmfReadGrid xdmfreadgrid_
55 #define XdmfReadGridAtIndex xdmfreadgridatindex_
56 #define XdmfGetNumberOfGrids xdmfgetnumberofgrids_
57 #define XdmfGetNumberOfPoints xdmfgetnumberofpoints_
58 #define XdmfReadPointValues xdmfreadpointvalues_
59 #define XdmfGetNumberOfAttributeValues xdmfgetnumberofattributevalues_
60 #define XdmfReadAttributeValues xdmfreadattributevalues_
61 #define XdmfReadInformationValue xdmfreadinformationvalue_
62 #define XdmfGetTime xdmfgettime_
63 #define XdmfWriteGrid xdmfwritegrid_
64 #define XdmfWriteToFile xdmfwritetofile_
65 #define XdmfSerialize xdmfserialize_
66 #define XdmfGetDOM xdmfgetdom_
67 #define XdmfClose xdmfclose_
68 #endif
69 /**
70  *
71  * Initialize a new Xdmf file.
72  *
73  */
XdmfFortran(char * outputName)74 XdmfFortran::XdmfFortran(char * outputName)
75 {
76   myDOM = new XdmfDOM();
77   myRoot = new XdmfRoot();
78   myDomain = new XdmfDomain();
79   myTopology = NULL;
80   myGeometry = NULL;
81   currentTime = -1;
82 
83   myRoot->SetDOM(myDOM);
84   myRoot->Build();
85   myRoot->Insert(myDomain);
86   myName = outputName;
87 }
88 
~XdmfFortran()89 XdmfFortran::~XdmfFortran()
90 {
91   this->Destroy();
92 }
93 
94 void
Destroy()95 XdmfFortran::Destroy()
96 {
97   currentTime = -1;
98 
99   delete myGeometry;
100   myGeometry = NULL;
101   delete myTopology;
102   myTopology = NULL;
103 
104   while (!myAttributes.empty())
105   {
106     delete myAttributes.back();
107     myAttributes.pop_back();
108   }
109 
110   while (!myInformations.empty())
111   {
112     delete myInformations.back();
113     myInformations.pop_back();
114   }
115 
116   while (!myCollections.empty())
117   {
118     delete myCollections.top();
119     myCollections.pop();
120   }
121 
122   delete myDOM;
123   delete myRoot;
124   delete myDomain;
125 }
126 
127 /**
128  *
129  * Set a time to be assigned to the next grid.
130  *
131  */
132 void
SetTime(double * t)133 XdmfFortran::SetTime(double * t)
134 {
135   currentTime = *t;
136 }
137 
138 /**
139  *
140  * Add a collection to the XdmfDOM.  Collections can be 'Spatial' or 'Temporal' type.
141  * Nested collections are supported.
142  *
143  */
144 void
AddCollection(char * collectionType)145 XdmfFortran::AddCollection(char * collectionType)
146 {
147   XdmfGrid * currentCollection = new XdmfGrid();
148   currentCollection->SetGridType(XDMF_GRID_COLLECTION);
149   currentCollection->SetCollectionTypeFromString(collectionType);
150   if (myCollections.empty())
151   {
152     myDomain->Insert(currentCollection);
153   }
154   else
155   {
156     myCollections.top()->Insert(currentCollection);
157   }
158   currentCollection->Build();
159   myCollections.push(currentCollection);
160 }
161 
162 /**
163  *
164  * Close the current open collection.  If within a nested collection, close
165  * the most deeply nested collection.
166  *
167  */
168 void
CloseCollection()169 XdmfFortran::CloseCollection()
170 {
171   if (!myCollections.empty())
172   {
173     delete myCollections.top();
174     myCollections.pop();
175   }
176 }
177 
178 /**
179  *
180  * Set the topology type to be assigned to the next grid.
181  * Only XDMF_INT_32 type currently supported for Topology --> INTEGER*4
182  *
183  */
184 void
SetGridTopology(char * topologyType,int * numberOfElements,XdmfInt32 * conns)185 XdmfFortran::SetGridTopology(char * topologyType, int * numberOfElements, XdmfInt32 * conns)
186 {
187   std::stringstream shape;
188   shape << *numberOfElements;
189   this->SetGridTopologyFromShape(topologyType, (char*) shape.str().c_str(), conns);
190 }
191 
192 /**
193  *
194  * Set the topology type to be assigned to the next grid.
195  * Only XDMF_INT_32 type currently supported for Topology --> INTEGER*4
196  *
197  */
198 void
SetGridTopologyFromShape(char * topologyType,char * shape,XdmfInt32 * conns)199 XdmfFortran::SetGridTopologyFromShape(char * topologyType, char * shape, XdmfInt32 * conns)
200 {
201   myTopology = new XdmfTopology();
202   myTopology->SetTopologyTypeFromString(topologyType);
203   myTopology->GetShapeDesc()->SetShapeFromString(shape);
204 
205   if (myTopology->GetClass() != XDMF_STRUCTURED && myTopology->GetTopologyType() != XDMF_POLYVERTEX)
206   {
207 
208     XdmfArray * myConnections = myTopology->GetConnectivity();
209     myConnections->SetNumberType(XDMF_INT32_TYPE);
210     myConnections->SetNumberOfElements(myTopology->GetNumberOfElements()
211         * myTopology->GetNodesPerElement());
212     myConnections->SetValues(0, conns, myTopology->GetNumberOfElements()
213         * myTopology->GetNodesPerElement());
214   }
215 }
216 
217 /**
218  *
219  * Set the geometry type to be assigned to the next grid.
220  *
221  */
222 void
SetGridGeometry(char * geometryType,char * numberType,int * numberOfPoints,XdmfPointer * points)223 XdmfFortran::SetGridGeometry(char * geometryType, char * numberType, int * numberOfPoints,
224     XdmfPointer * points)
225 {
226   myGeometry = new XdmfGeometry();
227   myGeometry->SetGeometryTypeFromString(geometryType);
228   myGeometry->SetNumberOfPoints(*numberOfPoints);
229 
230   XdmfArray * myPoints = myGeometry->GetPoints();
231   myPoints->SetNumberTypeFromString(numberType);
232 
233   switch (myGeometry->GetGeometryType())
234     {
235   case XDMF_GEOMETRY_XYZ:
236     myPoints->SetNumberOfElements(*numberOfPoints * 3);
237     break;
238   case XDMF_GEOMETRY_X_Y_Z:
239     myPoints->SetNumberOfElements(*numberOfPoints * 3);
240     break;
241   case XDMF_GEOMETRY_XY:
242     myPoints->SetNumberOfElements(*numberOfPoints * 2);
243     break;
244   case XDMF_GEOMETRY_X_Y:
245     myPoints->SetNumberOfElements(*numberOfPoints * 2);
246     break;
247   case XDMF_GEOMETRY_VXVYVZ:
248     //TODO: FIX THIS
249     myPoints->SetNumberOfElements(*numberOfPoints * 3);
250     break;
251   case XDMF_GEOMETRY_ORIGIN_DXDYDZ:
252     myPoints->SetNumberOfElements(6);
253     break;
254   default:
255     myPoints->SetNumberOfElements(*numberOfPoints * 3);
256     break;
257     }
258   WriteToXdmfArray(myPoints, points);
259 }
260 
261 /**
262  *
263  * Add an attribute to be written to the next grid.  Multiple attributes can
264  * be added and written to a single grid.
265  *
266  */
267 void
AddGridAttribute(char * attributeName,char * numberType,char * attributeCenter,char * attributeType,int * numberOfPoints,XdmfPointer * data)268 XdmfFortran::AddGridAttribute(char * attributeName, char * numberType, char * attributeCenter,
269     char * attributeType, int * numberOfPoints, XdmfPointer * data)
270 {
271   XdmfAttribute * currAttribute = new XdmfAttribute();
272   currAttribute->SetName(attributeName);
273   currAttribute->SetAttributeCenterFromString(attributeCenter);
274   currAttribute->SetAttributeTypeFromString(attributeType);
275   currAttribute->SetDeleteOnGridDelete(true);
276 
277   XdmfArray * array = currAttribute->GetValues();
278   array->SetNumberTypeFromString(numberType);
279   array->SetNumberOfElements(*numberOfPoints);
280   WriteToXdmfArray(array, data);
281   myAttributes.push_back(currAttribute);
282 }
283 
284 /**
285  *
286  * Add an attribute with shape to be written to the next grid.  Multiple attributes can
287  * be added and written to a single grid. Dimensions are specified using shape string rather than int.
288  * Extra argument supports setting units string of attribute.
289  *
290  */
291 void
AddGridAttributeFromShape(char * attributeName,char * numberType,char * attributeCenter,char * attributeType,char * shape,char * units,XdmfPointer * data)292 XdmfFortran::AddGridAttributeFromShape(char * attributeName, char * numberType,
293     char * attributeCenter, char * attributeType, char * shape, char * units, XdmfPointer * data)
294 {
295   XdmfAttribute * currAttribute = new XdmfAttribute();
296   currAttribute->SetName(attributeName);
297   currAttribute->SetUnits(units);
298   currAttribute->SetAttributeCenterFromString(attributeCenter);
299   currAttribute->SetAttributeTypeFromString(attributeType);
300   currAttribute->SetDeleteOnGridDelete(true);
301 
302   XdmfArray * array = currAttribute->GetValues();
303   array->SetNumberTypeFromString(numberType);
304   array->SetShapeFromString(shape); // dims via shape string
305   WriteToXdmfArray(array, data);
306   myAttributes.push_back(currAttribute);
307 }
308 
309 /**
310  *
311  * Add an attribute to be written to the next grid.  Multiple attributes can
312  * be added and written to a single grid.
313  *
314  */
315 void
AddGridInformation(char * informationName,char * value)316 XdmfFortran::AddGridInformation(char * informationName, char * value)
317 {
318   XdmfInformation * currInfo = new XdmfInformation();
319   currInfo->SetName(informationName);
320   currInfo->SetValue(value);
321   currInfo->SetDeleteOnGridDelete(true);
322 
323   myInformations.push_back(currInfo);
324 }
325 
326 /**
327  *
328  * Add an attribute to the current collection.  If we are not within a collection do nothing!
329  *
330  */
331 void
AddCollectionAttribute(char * attributeName,char * numberType,char * attributeCenter,char * attributeType,int * numberOfPoints,XdmfPointer * data)332 XdmfFortran::AddCollectionAttribute(char * attributeName, char * numberType,
333     char * attributeCenter, char * attributeType, int * numberOfPoints, XdmfPointer * data)
334 {
335   if (!myCollections.empty())
336   {
337     XdmfAttribute * currAttribute = new XdmfAttribute();
338     currAttribute->SetName(attributeName);
339     currAttribute->SetAttributeCenterFromString(attributeCenter);
340     currAttribute->SetAttributeTypeFromString(attributeType);
341     currAttribute->SetDeleteOnGridDelete(true);
342 
343     XdmfArray * array = currAttribute->GetValues();
344     array->SetNumberTypeFromString(numberType);
345     array->SetNumberOfElements(*numberOfPoints);
346     WriteToXdmfArray(array, data);
347     myCollections.top()->Insert(currAttribute);
348     myCollections.top()->Build();
349   }
350 }
351 
352 /**
353  *
354  * Add an attribute to the current collection.  If we are not within a collection add to the top level domain
355  *
356  */
357 void
AddCollectionInformation(char * informationName,char * value)358 XdmfFortran::AddCollectionInformation(char * informationName, char * value)
359 {
360   XdmfInformation * currInfo = new XdmfInformation();
361   currInfo->SetName(informationName);
362   currInfo->SetValue(value);
363   currInfo->SetDeleteOnGridDelete(true);
364   if (!myCollections.empty())
365   {
366 
367     myCollections.top()->Insert(currInfo);
368     myCollections.top()->Build();
369   }
370   else
371   {
372     myDomain->Insert(currInfo);
373     myDomain->Build();
374   }
375 }
376 
377 /**
378  *
379  * Write out "generic" data to XDMF.  This writes out data to the end of the top-level domain or the current collection.  It is independent of any grids.
380  * Currently supports only writing a single dataitem.
381  *
382  */
383 void
AddArray(char * name,char * numberType,int * numberOfValues,XdmfPointer * data)384 XdmfFortran::AddArray(char * name, char * numberType, int * numberOfValues, XdmfPointer * data)
385 {
386   XdmfSet * currSet = new XdmfSet();
387   currSet->SetDOM(myDOM);
388   currSet->SetSetType(XDMF_SET_TYPE_NODE);
389   currSet->SetName(name);
390   currSet->SetDeleteOnGridDelete(true);
391 
392   // Copy Elements from Set to XdmfArray
393   XdmfArray * array = currSet->GetIds();
394   array->SetNumberTypeFromString(numberType);
395   array->SetNumberOfElements(*numberOfValues);
396   std::stringstream heavyDataName;
397   heavyDataName << myName << ".h5:/" << name;
398   ;
399   array->SetHeavyDataSetName(heavyDataName.str().c_str());
400   WriteToXdmfArray(array, data);
401 
402   if (!myCollections.empty())
403   {
404     myCollections.top()->Insert(currSet);
405     currSet->Build();
406   }
407   else
408   {
409     XdmfGrid * myGrid = new XdmfGrid();
410     myGrid->SetDOM(myDOM);
411     myGrid->SetElement(myDOM->FindElement("Domain"));
412     myGrid->Insert(currSet);
413     currSet->Build();
414     delete myGrid;
415   }
416 }
417 
418 /**
419  *
420  * Read an Xdmf file into the current XdmfDOM.  Must call XdmfReadGrid() to read in associated geometry
421  * topology, and attributes.
422  *
423  */
424 void
ReadFile(char * filePath)425 XdmfFortran::ReadFile(char * filePath)
426 {
427   // Clear state and start over before reading file
428   this->Destroy();
429 
430   myDOM = new XdmfDOM();
431   myRoot = new XdmfRoot();
432   myDomain = new XdmfDomain();
433 
434   myDOM->Parse(filePath);
435   myDomain->SetElement(myDOM->FindElement("Domain"));
436   myRoot->SetElement(myDOM->GetRoot());
437 
438   // Perhaps we should support collections more on this part?
439   while (!myCollections.empty())
440   {
441     delete myCollections.top();
442     myCollections.pop();
443   }
444   myGridPaths.clear();
445   myGridNames.clear();
446   this->ReadFilePriv(myDomain->GetElement());
447 }
448 
449 void
ReadFilePriv(XdmfXmlNode currElement)450 XdmfFortran::ReadFilePriv(XdmfXmlNode currElement)
451 {
452   XdmfGrid currGrid = XdmfGrid();
453   for (int i = 0; i < myDOM->FindNumberOfElements("Grid", currElement); i++)
454   {
455     currGrid.SetDOM(myDOM);
456     currGrid.SetElement(myDOM->FindElement("Grid", i, currElement));
457     currGrid.Update();
458     if (currGrid.GetGridType() != XDMF_GRID_COLLECTION)
459     {
460       myGridPaths.push_back(myDOM->GetPath(currGrid.GetElement()));
461       std::string gridName = currGrid.GetName();
462       if (gridName.find_last_of("_") != std::string::npos)
463       {
464         try
465         {
466           atoi(gridName.substr(gridName.find_last_of("_") + 1, gridName.size()
467               - gridName.find_last_of("_")).c_str());
468           gridName = gridName.substr(0, gridName.find_last_of("_"));
469         }
470         catch (int)
471         {
472         }
473       }
474       if (myGridNames.find(gridName.c_str()) == myGridNames.end())
475       {
476         myGridNames[gridName.c_str()] = 1;
477       }
478       else
479       {
480         myGridNames[gridName.c_str()]++;
481       }
482     }
483     this->ReadFilePriv(currGrid.GetElement());
484   }
485 }
486 
487 /**
488  *
489  * Read a grid in the current XdmfDOM into XdmfGeometry, XdmfTopology, and XdmfAttribute elements.
490  * An XdmfReadGrid() followed by a XdmfWriteGrid() will make a copy of the grid.
491  *
492  */
493 void
ReadGrid(char * gridName)494 XdmfFortran::ReadGrid(char * gridName)
495 {
496   ReadGridPriv(gridName, myDomain->GetElement());
497 }
498 
499 /**
500  *
501  * Read a grid in the current XdmfDOM into XdmfGeometry, XdmfTopology, and XdmfAttribute elements.
502  * An XdmfReadGrid() followed by a XdmfWriteGrid() will make a copy of the grid.
503  *
504  */
505 void
ReadGridAtIndex(int * gridIndex)506 XdmfFortran::ReadGridAtIndex(int * gridIndex)
507 {
508   ReadGridPriv(myGridPaths[*gridIndex].c_str());
509 }
510 
511 /**
512  *
513  * Helper function for XdmfReadGrid.  Ensures that all grids are traversed and that the method works
514  * even within collections.
515  *
516  */
517 void
ReadGridPriv(char * gridName,XdmfXmlNode currElement)518 XdmfFortran::ReadGridPriv(char * gridName, XdmfXmlNode currElement)
519 {
520   XdmfGrid currGrid = XdmfGrid();
521   for (int i = 0; i < myDOM->FindNumberOfElements("Grid", currElement); i++)
522   {
523     currGrid.SetDOM(myDOM);
524     currGrid.SetElement(myDOM->FindElement("Grid", i, currElement));
525     currGrid.Update();
526     if (currGrid.GetGridType() != XDMF_GRID_COLLECTION)
527     {
528       if (strcmp(gridName, currGrid.GetName()) == 0)
529       {
530         return this->ReadGridPriv(myDOM->GetPath(currGrid.GetElement()));
531       }
532     }
533     this->ReadGridPriv(gridName, currGrid.GetElement());
534   }
535 }
536 
537 /**
538  *
539  * Helper function for XdmfReadGrid.  Ensures that all grids are traversed and that the method works
540  * even within collections.
541  *
542  */
543 void
ReadGridPriv(XdmfConstString gridPath)544 XdmfFortran::ReadGridPriv(XdmfConstString gridPath)
545 {
546   XdmfGrid currGrid = XdmfGrid();
547   currGrid.SetDOM(myDOM);
548   currGrid.SetElement(myDOM->FindElementByPath(gridPath));
549   currGrid.Update();
550 
551   delete myGeometry;
552   delete myTopology;
553 
554   myGeometry = new XdmfGeometry();
555   myGeometry->SetGeometryType(currGrid.GetGeometry()->GetGeometryType());
556   myGeometry->SetNumberOfPoints(currGrid.GetGeometry()->GetNumberOfPoints());
557   myGeometry->SetPoints(currGrid.GetGeometry()->GetPoints()->Clone());
558 
559   myTopology = new XdmfTopology();
560   myTopology->SetTopologyType(currGrid.GetTopology()->GetTopologyType());
561   myTopology->SetNumberOfElements(currGrid.GetTopology()->GetNumberOfElements());
562   myTopology->SetConnectivity(currGrid.GetTopology()->GetConnectivity()->Clone());
563 
564   while (!myAttributes.empty())
565   {
566     delete myAttributes.back();
567     myAttributes.pop_back();
568   }
569 
570   for (int j = 0; j < currGrid.GetNumberOfAttributes(); j++)
571   {
572     currGrid.GetAttribute(j)->Update();
573     XdmfAttribute * currAttribute = new XdmfAttribute();
574     currAttribute->SetName(currGrid.GetAttribute(j)->GetName());
575     currAttribute->SetAttributeCenter(currGrid.GetAttribute(j)->GetAttributeCenter());
576     currAttribute->SetAttributeType(currGrid.GetAttribute(j)->GetAttributeType());
577     currAttribute->SetDeleteOnGridDelete(true);
578     currAttribute->SetValues(currGrid.GetAttribute(j)->GetValues()->Clone());
579     myAttributes.push_back(currAttribute);
580   }
581 
582   for (int j = 0; j < currGrid.GetNumberOfInformations(); j++)
583   {
584     currGrid.GetInformation(j)->UpdateInformation();
585     XdmfInformation * currInformation = new XdmfInformation();
586     currInformation->SetName(currGrid.GetInformation(j)->GetName());
587     currInformation->SetValue(currGrid.GetInformation(j)->GetValue());
588     currInformation->SetDeleteOnGridDelete(true);
589     myInformations.push_back(currInformation);
590   }
591 }
592 
593 /**
594  *
595  * Returns the number of grids in the current open file.  This ignores collections.
596  *
597  */
598 void
GetNumberOfGrids(XdmfInt32 * toReturn)599 XdmfFortran::GetNumberOfGrids(XdmfInt32 * toReturn)
600 {
601   *toReturn = myGridPaths.size();
602 }
603 
604 /**
605  *
606  * Returns the number of points in the current open grid (the current active XdmfGeometry Element).  This is either
607  * from a current read-in file or from a created but unwritten grid.  If no geometry element is present, return -1.
608  *
609  */
610 void
GetNumberOfPoints(XdmfInt32 * toReturn)611 XdmfFortran::GetNumberOfPoints(XdmfInt32 * toReturn)
612 {
613   if (myGeometry != NULL)
614   {
615     *toReturn = (XdmfInt32) myGeometry->GetNumberOfPoints();
616   }
617   else
618   {
619     *toReturn = -1;
620   }
621 }
622 
623 /**
624  *
625  * Reads the point values from the current geometry into the passed array pointer.  If the geometry has not been created
626  * no values are read.
627  *
628  */
629 void
ReadPointValues(char * numberType,XdmfInt32 * startIndex,XdmfPointer * arrayToFill,XdmfInt32 * numberOfValues,XdmfInt32 * arrayStride,XdmfInt32 * valuesStride)630 XdmfFortran::ReadPointValues(char * numberType, XdmfInt32 * startIndex, XdmfPointer * arrayToFill,
631     XdmfInt32 * numberOfValues, XdmfInt32 * arrayStride, XdmfInt32 * valuesStride)
632 {
633   if (myGeometry != NULL)
634   {
635     this->ReadFromXdmfArray(myGeometry->GetPoints(), numberType, startIndex, arrayToFill,
636         numberOfValues, arrayStride, valuesStride);
637   }
638 }
639 
640 /**
641  *
642  * Returns the number of values in the specified attribute.  Iterates over all current open attributes to find the
643  * specified attribute name and returns the number of values it contains.  If no attribute is found, return -1.
644  *
645  */
646 void
GetNumberOfAttributeValues(char * attributeName,XdmfInt32 * toReturn)647 XdmfFortran::GetNumberOfAttributeValues(char * attributeName, XdmfInt32 * toReturn)
648 {
649   for (unsigned int i = 0; i < myAttributes.size(); i++)
650   {
651     if (strcmp(myAttributes[i]->GetName(), attributeName) == 0)
652     {
653       *toReturn = (XdmfInt32) myAttributes[i]->GetValues()->GetNumberOfElements();
654       return;
655     }
656   }
657   *toReturn = -1;
658 }
659 
660 /**
661  *
662  * Reads the values from the specified attribute into the passed array pointer.  If the attribute cannot be found,
663  * no values are read.
664  *
665  */
666 void
ReadAttributeValues(char * attributeName,char * numberType,XdmfInt32 * startIndex,XdmfPointer * arrayToFill,XdmfInt32 * numberOfValues,XdmfInt32 * arrayStride,XdmfInt32 * valuesStride)667 XdmfFortran::ReadAttributeValues(char * attributeName, char * numberType, XdmfInt32 * startIndex,
668     XdmfPointer * arrayToFill, XdmfInt32 * numberOfValues, XdmfInt32 * arrayStride,
669     XdmfInt32 * valuesStride)
670 {
671   for (unsigned int i = 0; i < myAttributes.size(); i++)
672   {
673     if (strcmp(myAttributes[i]->GetName(), attributeName) == 0)
674     {
675       this->ReadFromXdmfArray(myAttributes[i]->GetValues(), numberType, startIndex, arrayToFill,
676           numberOfValues, arrayStride, valuesStride);
677     }
678   }
679 }
680 
681 /**
682  *
683  * Reads the values from the specified information element into the passed const string pointer.  If the information element cannot be found,
684  * no values are passed.  Information elements at the top level domain are searched first, followed by the currently loaded grid.
685  *
686  */
687 void
ReadInformationValue(char * informationName,char * toReturn)688 XdmfFortran::ReadInformationValue(char * informationName, char * toReturn)
689 {
690   // TODO: Make this work better for collections as well!
691   for (unsigned int i = 0; i < myInformations.size(); i++)
692   {
693     if (strcmp(informationName, myInformations[i]->GetName()) == 0)
694     {
695       strcpy(toReturn, myInformations[i]->GetValue());
696       return;
697     }
698   }
699 
700   for (int i = 0; i < myDOM->FindNumberOfElements("Information", myDomain->GetElement()); i++)
701   {
702     XdmfInformation currInfo = XdmfInformation();
703     currInfo.SetDOM(myDOM);
704     currInfo.SetElement(myDOM->FindElement("Information", i, myDomain->GetElement(), 0));
705     currInfo.UpdateInformation();
706     if (strcmp(informationName, currInfo.GetName()) == 0)
707     {
708       strcpy(toReturn, currInfo.GetValue());
709       return;
710     }
711   }
712 }
713 
714 /**
715  *
716  * Return the currentTime
717  *
718  */
719 void
GetTime(XdmfFloat64 * toReturn)720 XdmfFortran::GetTime(XdmfFloat64 * toReturn)
721 {
722   *toReturn = currentTime;
723 }
724 
725 /**
726  *
727  * Add a grid to the XdmfDOM.  Assign the current topology, geometry, and grid attributes
728  * to grid.  If within a collection, add grid to the collection, otherwise
729  * add to the top level domain.  Assign time value if value is nonnegative.
730  *
731  */
732 void
WriteGrid(char * gridName)733 XdmfFortran::WriteGrid(char * gridName)
734 {
735   XdmfGrid * grid = new XdmfGrid();
736 
737   if (myTopology == NULL)
738   {
739     cout << "Must set a topology before the grid can be written" << endl;
740     delete grid;
741     return;
742   }
743 
744   if (myGeometry == NULL)
745   {
746     cout << "Must set a geometry before the grid can be written" << endl;
747     delete grid;
748     return;
749   }
750 
751   // If we try to write over the same grid, modify the grid name...
752   std::stringstream totalGridName;
753   if (myGridNames.find(gridName) == myGridNames.end())
754   {
755     myGridNames[gridName] = 1;
756     totalGridName << gridName;
757   }
758   else
759   {
760     myGridNames[gridName]++;
761     totalGridName << gridName << "_" << myGridNames[gridName];
762   }
763 
764   grid->SetName(totalGridName.str().c_str());
765 
766   // Set Topology
767   if (myTopology->GetClass() != XDMF_STRUCTURED && myTopology->GetTopologyType() != XDMF_POLYVERTEX)
768   {
769     //Modify HDF5 names so we aren't writing over top of our data!
770     std::stringstream topologyDataName;
771     topologyDataName << myName << ".h5:/" << totalGridName.str() << "/Connections";
772     myTopology->GetConnectivity()->SetHeavyDataSetName(topologyDataName.str().c_str());
773   }
774   grid->SetTopology(myTopology);
775 
776   // Set Geometry
777   std::stringstream geometryDataName;
778   geometryDataName << myName << ".h5:/" << totalGridName.str() << "/XYZ";
779   myGeometry->GetPoints()->SetHeavyDataSetName(geometryDataName.str().c_str());
780   grid->SetGeometry(myGeometry);
781 
782   if (!myCollections.empty())
783   {
784     myCollections.top()->Insert(grid);
785   }
786   else
787   {
788     myDomain->Insert(grid);
789   }
790 
791   XdmfTime * t = new XdmfTime();
792   if (currentTime >= 0)
793   {
794     t->SetTimeType(XDMF_TIME_SINGLE);
795     t->SetValue(currentTime);
796     grid->Insert(t);
797     currentTime = -1;
798   }
799 
800   while (myInformations.size() > 0)
801   {
802     grid->Insert(myInformations.back());
803     myInformations.pop_back();
804   }
805 
806   while (myAttributes.size() > 0)
807   {
808     XdmfAttribute * currAttribute = myAttributes.back();
809     std::stringstream attributeDataName;
810     attributeDataName << myName << ".h5:/" << totalGridName.str() << "/"
811         << currAttribute->GetName();
812     currAttribute->GetValues()->SetHeavyDataSetName(attributeDataName.str().c_str());
813 
814     grid->Insert(currAttribute);
815     myAttributes.pop_back();
816   }
817 
818   grid->Build();
819 
820   myGridPaths.push_back(myDOM->GetPath(grid->GetElement()));
821 
822   // If we are within a collection this will be deleted on collection deletion
823   if (myCollections.empty())
824   {
825     delete grid;
826   }
827   myTopology = NULL;
828   myGeometry = NULL;
829 }
830 
831 /**
832  *
833  * Write constructed Xdmf file to disk with filename created upon initialization
834  *
835  */
836 void
WriteToFile()837 XdmfFortran::WriteToFile()
838 {
839   std::stringstream dataName;
840   dataName << myName << ".xmf";
841   myDOM->Write(dataName.str().c_str());
842 }
843 
844 /**
845  *
846  * Print current XdmfDOM to console
847  *
848  */
849 void
Serialize()850 XdmfFortran::Serialize()
851 {
852   cout << myDOM->Serialize() << endl;
853 }
854 
855 /**
856  *
857  * Copy current XdmfDOM to memory pointed to by charPointer
858  *
859  */
860 void
GetDOM(char * charPointer)861 XdmfFortran::GetDOM(char * charPointer)
862 {
863   strcpy(charPointer, myDOM->Serialize());
864 }
865 
866 /**
867  *
868  * Helper function to write different datatypes to an XdmfArray.
869  *
870  */
871 void
WriteToXdmfArray(XdmfArray * array,XdmfPointer * data)872 XdmfFortran::WriteToXdmfArray(XdmfArray * array, XdmfPointer * data)
873 {
874   switch (array->GetNumberType())
875     {
876   case XDMF_INT8_TYPE:
877     array->SetValues(0, (XdmfInt8*) data, array->GetNumberOfElements());
878     return;
879   case XDMF_INT16_TYPE:
880     array->SetValues(0, (XdmfInt16*) data, array->GetNumberOfElements());
881     return;
882   case XDMF_INT32_TYPE:
883     array->SetValues(0, (XdmfInt32*) data, array->GetNumberOfElements());
884     return;
885   case XDMF_INT64_TYPE:
886     array->SetValues(0, (XdmfInt64*) data, array->GetNumberOfElements());
887     return;
888   case XDMF_FLOAT32_TYPE:
889     array->SetValues(0, (XdmfFloat32*) data, array->GetNumberOfElements());
890     return;
891   case XDMF_FLOAT64_TYPE:
892     array->SetValues(0, (XdmfFloat64*) data, array->GetNumberOfElements());
893     return;
894   case XDMF_UINT8_TYPE:
895     array->SetValues(0, (XdmfUInt8*) data, array->GetNumberOfElements());
896     return;
897   case XDMF_UINT16_TYPE:
898     array->SetValues(0, (XdmfUInt16*) data, array->GetNumberOfElements());
899     return;
900   case XDMF_UINT32_TYPE:
901     array->SetValues(0, (XdmfUInt32*) data, array->GetNumberOfElements());
902     return;
903   default:
904     array->SetValues(0, (XdmfFloat64*) data, array->GetNumberOfElements());
905     return;
906     }
907 }
908 
909 /**
910  *
911  * Helper function to read different datatypes from an XdmfArray.
912  *
913  */
914 void
ReadFromXdmfArray(XdmfArray * array,char * numberType,XdmfInt32 * startIndex,XdmfPointer * arrayToFill,XdmfInt32 * numberOfValues,XdmfInt32 * arrayStride,XdmfInt32 * valuesStride)915 XdmfFortran::ReadFromXdmfArray(XdmfArray * array, char * numberType, XdmfInt32 * startIndex,
916     XdmfPointer * arrayToFill, XdmfInt32 * numberOfValues, XdmfInt32 * arrayStride,
917     XdmfInt32 * valuesStride)
918 {
919   XdmfArray a = XdmfArray();
920   a.SetNumberTypeFromString(numberType);
921   switch (a.GetNumberType())
922     {
923   case XDMF_INT8_TYPE:
924     array->GetValues(*startIndex, (XdmfInt8*) arrayToFill, *numberOfValues, *arrayStride,
925         *valuesStride);
926     return;
927   case XDMF_INT16_TYPE:
928     array->GetValues(*startIndex, (XdmfInt16*) arrayToFill, *numberOfValues, *arrayStride,
929         *valuesStride);
930     return;
931   case XDMF_INT32_TYPE:
932     array->GetValues(*startIndex, (XdmfInt32*) arrayToFill, *numberOfValues, *arrayStride,
933         *valuesStride);
934     return;
935   case XDMF_INT64_TYPE:
936     array->GetValues(*startIndex, (XdmfInt64*) arrayToFill, *numberOfValues, *arrayStride,
937         *valuesStride);
938     return;
939   case XDMF_FLOAT32_TYPE:
940     array->GetValues(*startIndex, (XdmfFloat32*) arrayToFill, *numberOfValues, *arrayStride,
941         *valuesStride);
942     return;
943   case XDMF_FLOAT64_TYPE:
944     array->GetValues(*startIndex, (XdmfFloat64*) arrayToFill, *numberOfValues, *arrayStride,
945         *valuesStride);
946     return;
947   case XDMF_UINT8_TYPE:
948     array->GetValues(*startIndex, (XdmfUInt8*) arrayToFill, *numberOfValues, *arrayStride,
949         *valuesStride);
950     return;
951   case XDMF_UINT16_TYPE:
952     array->GetValues(*startIndex, (XdmfUInt16*) arrayToFill, *numberOfValues, *arrayStride,
953         *valuesStride);
954     return;
955   case XDMF_UINT32_TYPE:
956     array->GetValues(*startIndex, (XdmfUInt32*) arrayToFill, *numberOfValues, *arrayStride,
957         *valuesStride);
958     return;
959   default:
960     array->GetValues(*startIndex, (XdmfFloat64*) arrayToFill, *numberOfValues, *arrayStride,
961         *valuesStride);
962     return;
963     }
964 }
965 
966 //
967 // C++ will mangle the name based on the argument list. This tells the
968 // compiler not to mangle the name so we can call it from 'C' (but
969 // really Fortran in this case)
970 //
971 extern "C"
972 {
973 
974   /**
975    *
976    * Initialize a new Xdmf file.
977    *
978    */
979   void XDMF_UTILS_DLL
XdmfInit(long * pointer,char * outputName)980   XdmfInit(long * pointer, char * outputName)
981   {
982     printf("Inside XdmfInit - outputName is: %s\n", outputName);
983 
984     XdmfFortran * myPointer = new XdmfFortran(outputName);
985     *pointer = (long) myPointer;
986   }
987 
988   void XDMF_UTILS_DLL
XdmfSetTime(long * pointer,double * t)989   XdmfSetTime(long * pointer, double * t)
990   {
991     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
992     myPointer->SetTime(t);
993   }
994 
995   void XDMF_UTILS_DLL
XdmfAddCollection(long * pointer,char * collectionType)996   XdmfAddCollection(long * pointer, char * collectionType)
997   {
998     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
999     myPointer->AddCollection(collectionType);
1000   }
1001 
1002   void XDMF_UTILS_DLL
XdmfCloseCollection(long * pointer)1003   XdmfCloseCollection(long * pointer)
1004   {
1005     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1006     myPointer->CloseCollection();
1007   }
1008 
1009   void XDMF_UTILS_DLL
XdmfSetGridTopology(long * pointer,char * topologyType,int * numberOfElements,XdmfInt32 * conns)1010   XdmfSetGridTopology(long * pointer, char * topologyType, int * numberOfElements,
1011       XdmfInt32 * conns)
1012   {
1013     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1014     myPointer->SetGridTopology(topologyType, numberOfElements, conns);
1015   }
1016 
1017   void XDMF_UTILS_DLL
XdmfSetGridTopologyFromShape(long * pointer,char * topologyType,char * shape,XdmfInt32 * conns)1018   XdmfSetGridTopologyFromShape(long * pointer, char * topologyType, char * shape, XdmfInt32 * conns)
1019   {
1020     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1021     myPointer->SetGridTopologyFromShape(topologyType, shape, conns);
1022   }
1023 
1024   void XDMF_UTILS_DLL
XdmfSetGridGeometry(long * pointer,char * geometryType,char * numberType,int * numberOfPoints,XdmfPointer * points)1025   XdmfSetGridGeometry(long * pointer, char * geometryType, char * numberType, int * numberOfPoints,
1026       XdmfPointer * points)
1027   {
1028     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1029     myPointer->SetGridGeometry(geometryType, numberType, numberOfPoints, points);
1030   }
1031 
1032   void XDMF_UTILS_DLL
XdmfAddGridAttribute(long * pointer,char * attributeName,char * numberType,char * attributeCenter,char * attributeType,int * numberOfPoints,XdmfPointer * data)1033   XdmfAddGridAttribute(long * pointer, char * attributeName, char * numberType,
1034       char * attributeCenter, char * attributeType, int * numberOfPoints, XdmfPointer * data)
1035   {
1036     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1037     myPointer->AddGridAttribute(attributeName, numberType, attributeCenter, attributeType,
1038         numberOfPoints, data);
1039   }
1040 
1041   void XDMF_UTILS_DLL
XdmfAddGridAttributeFromShape(long * pointer,char * attributeName,char * numberType,char * attributeCenter,char * attributeType,char * shape,char * units,XdmfPointer * data)1042   XdmfAddGridAttributeFromShape(long * pointer, char * attributeName, char * numberType,
1043       char * attributeCenter, char * attributeType, char * shape, char * units, XdmfPointer * data)
1044   {
1045     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1046     myPointer->AddGridAttributeFromShape(attributeName, numberType, attributeCenter, attributeType,
1047         shape, units, data);
1048   }
1049 
1050   void XDMF_UTILS_DLL
XdmfAddCollectionAttribute(long * pointer,char * attributeName,char * numberType,char * attributeCenter,char * attributeType,int * numberOfPoints,XdmfPointer * data)1051   XdmfAddCollectionAttribute(long * pointer, char * attributeName, char * numberType,
1052       char * attributeCenter, char * attributeType, int * numberOfPoints, XdmfPointer * data)
1053   {
1054     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1055     myPointer->AddCollectionAttribute(attributeName, numberType, attributeCenter, attributeType,
1056         numberOfPoints, data);
1057   }
1058 
1059   void XDMF_UTILS_DLL
XdmfAddGridInformation(long * pointer,char * informationName,char * value)1060   XdmfAddGridInformation(long * pointer, char * informationName, char * value)
1061   {
1062     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1063     myPointer->AddGridInformation(informationName, value);
1064   }
1065 
1066   void XDMF_UTILS_DLL
XdmfAddCollectionInformation(long * pointer,char * informationName,char * value)1067   XdmfAddCollectionInformation(long * pointer, char * informationName, char * value)
1068   {
1069     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1070     myPointer->AddCollectionInformation(informationName, value);
1071   }
1072 
1073   void XDMF_UTILS_DLL
XdmfAddArray(long * pointer,char * name,char * numberType,int * numberOfValues,XdmfPointer * data)1074   XdmfAddArray(long * pointer, char * name, char * numberType, int * numberOfValues,
1075       XdmfPointer * data)
1076   {
1077     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1078     myPointer->AddArray(name, numberType, numberOfValues, data);
1079   }
1080 
1081   void XDMF_UTILS_DLL
XdmfReadFile(long * pointer,char * filePath)1082   XdmfReadFile(long * pointer, char * filePath)
1083   {
1084     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1085     myPointer->ReadFile(filePath);
1086   }
1087 
1088   void XDMF_UTILS_DLL
XdmfReadGrid(long * pointer,char * gridName)1089   XdmfReadGrid(long * pointer, char * gridName)
1090   {
1091     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1092     myPointer->ReadGrid(gridName);
1093   }
1094 
1095   void XDMF_UTILS_DLL
XdmfReadGridAtIndex(long * pointer,int * gridIndex)1096   XdmfReadGridAtIndex(long * pointer, int * gridIndex)
1097   {
1098     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1099     myPointer->ReadGridAtIndex(gridIndex);
1100   }
1101 
1102   void XDMF_UTILS_DLL
XdmfGetNumberOfGrids(long * pointer,XdmfInt32 * toReturn)1103   XdmfGetNumberOfGrids(long * pointer, XdmfInt32 * toReturn)
1104   {
1105     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1106     myPointer->GetNumberOfGrids(toReturn);
1107   }
1108 
1109   void XDMF_UTILS_DLL
XdmfGetNumberOfPoints(long * pointer,XdmfInt32 * toReturn)1110   XdmfGetNumberOfPoints(long * pointer, XdmfInt32 * toReturn)
1111   {
1112     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1113     myPointer->GetNumberOfPoints(toReturn);
1114   }
1115 
1116   void XDMF_UTILS_DLL
XdmfReadPointValues(long * pointer,char * numberType,XdmfInt32 * startIndex,XdmfPointer * arrayToFill,XdmfInt32 * numberOfValues,XdmfInt32 * arrayStride,XdmfInt32 * valuesStride)1117   XdmfReadPointValues(long * pointer, char * numberType, XdmfInt32 * startIndex,
1118       XdmfPointer * arrayToFill, XdmfInt32 * numberOfValues, XdmfInt32 * arrayStride,
1119       XdmfInt32 * valuesStride)
1120   {
1121     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1122     myPointer->ReadPointValues(numberType, startIndex, arrayToFill, numberOfValues, arrayStride,
1123         valuesStride);
1124   }
1125 
1126   void XDMF_UTILS_DLL
XdmfGetNumberOfAttributeValues(long * pointer,char * attributeName,XdmfInt32 * toReturn)1127   XdmfGetNumberOfAttributeValues(long * pointer, char * attributeName, XdmfInt32 * toReturn)
1128   {
1129     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1130     myPointer->GetNumberOfAttributeValues(attributeName, toReturn);
1131   }
1132 
1133   void XDMF_UTILS_DLL
XdmfReadAttributeValues(long * pointer,char * attributeName,char * numberType,XdmfInt32 * startIndex,XdmfPointer * arrayToFill,XdmfInt32 * numberOfValues,XdmfInt32 * arrayStride,XdmfInt32 * valuesStride)1134   XdmfReadAttributeValues(long * pointer, char * attributeName, char * numberType,
1135       XdmfInt32 * startIndex, XdmfPointer * arrayToFill, XdmfInt32 * numberOfValues,
1136       XdmfInt32 * arrayStride, XdmfInt32 * valuesStride)
1137   {
1138     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1139     myPointer->ReadAttributeValues(attributeName, numberType, startIndex, arrayToFill,
1140         numberOfValues, arrayStride, valuesStride);
1141   }
1142 
1143   void XDMF_UTILS_DLL
XdmfReadInformationValue(long * pointer,char * informationName,char * valueToReturn)1144   XdmfReadInformationValue(long * pointer, char * informationName, char * valueToReturn)
1145   {
1146     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1147     myPointer->ReadInformationValue(informationName, valueToReturn);
1148   }
1149 
1150   void XDMF_UTILS_DLL
XdmfGetTime(long * pointer,XdmfFloat64 * toReturn)1151   XdmfGetTime(long * pointer, XdmfFloat64 * toReturn)
1152   {
1153     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1154     myPointer->GetTime(toReturn);
1155   }
1156 
1157   void XDMF_UTILS_DLL
XdmfWriteGrid(long * pointer,char * gridName)1158   XdmfWriteGrid(long * pointer, char * gridName)
1159   {
1160     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1161     myPointer->WriteGrid(gridName);
1162   }
1163 
1164   void XDMF_UTILS_DLL
XdmfWriteToFile(long * pointer)1165   XdmfWriteToFile(long * pointer)
1166   {
1167     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1168     myPointer->WriteToFile();
1169   }
1170 
1171   void XDMF_UTILS_DLL
XdmfSerialize(long * pointer)1172   XdmfSerialize(long * pointer)
1173   {
1174     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1175     myPointer->Serialize();
1176   }
1177 
1178   void XDMF_UTILS_DLL
XdmfGetDOM(long * pointer,char * charPointer)1179   XdmfGetDOM(long * pointer, char * charPointer)
1180   {
1181     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1182     myPointer->GetDOM(charPointer);
1183   }
1184 
1185   /**
1186    *
1187    * Close XdmfFortran interface and clean memory
1188    *
1189    */
1190   void XDMF_UTILS_DLL
XdmfClose(long * pointer)1191   XdmfClose(long * pointer)
1192   {
1193     XdmfFortran * myPointer = (XdmfFortran *) *pointer;
1194     delete myPointer;
1195   }
1196 }
1197