VTK
vtkXMLWriterF.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkXMLWriterF.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
15 /*
16  * vtkXMLWriterF.h helps fortran programs call the C interface for
17  * writing VTK XML files. A program can use this by writing one
18  * vtkXMLWriterF.c file that includes this header. DO NOT INCLUDE
19  * THIS HEADER ELSEWHERE. The fortran program then compiles
20  * vtkXMLWriterF.c using a C compiler and links to the resulting
21  * object file.
22  */
23 
24 #if defined(__cplusplus)
25 # error "This should be included only by a .c file."
26 #endif
27 
28 /* Calls will be forwarded to the C interface. */
29 #include "vtkXMLWriterC.h"
30 
31 #include <stdio.h> /* fprintf */
32 #include <stdlib.h> /* malloc, free */
33 #include <string.h> /* memcpy */
34 
35 /* Define a static-storage default-zero-initialized table to store
36  writer objects for the fortran program. */
37 #define VTK_XMLWRITERF_MAX 256
39 
40 /* Fortran compilers expect certain symbol names for their calls to C
41  code. These macros help build the C symbols so that the fortran
42  program can link to them properly. The definitions here are
43  reasonable defaults but the source file that includes this can
44  define them appropriately for a particular compiler and override
45  these. */
46 #if !defined(VTK_FORTRAN_NAME)
47 # define VTK_FORTRAN_NAME(name, NAME) name##__
48 #endif
49 #if !defined(VTK_FORTRAN_ARG_STRING_POINTER)
50 # define VTK_FORTRAN_ARG_STRING_POINTER(name) const char* name##_ptr_arg
51 #endif
52 #if !defined(VTK_FORTRAN_ARG_STRING_LENGTH)
53 # define VTK_FORTRAN_ARG_STRING_LENGTH(name) , const long int name##_len_arg
54 #endif
55 #if !defined(VTK_FORTRAN_REF_STRING_POINTER)
56 # define VTK_FORTRAN_REF_STRING_POINTER(name) name##_ptr_arg
57 #endif
58 #if !defined(VTK_FORTRAN_REF_STRING_LENGTH)
59 # define VTK_FORTRAN_REF_STRING_LENGTH(name) ((int)name##_len_arg)
60 #endif
61 
62 /*--------------------------------------------------------------------------*/
63 /* vtkXMLWriterF_New */
64 void VTK_FORTRAN_NAME(vtkxmlwriterf_new, VTKXMLWRITERF_NEW)(
65  int* self
66  )
67 {
68  int i;
69 
70  /* Initialize result to failure. */
71  *self = 0;
72 
73  /* Search for a table entry to use for this object. */
74  for(i=1;i <= VTK_XMLWRITERF_MAX; ++i)
75  {
76  if(!vtkXMLWriterF_Table[i])
77  {
79  if(vtkXMLWriterF_Table[i])
80  {
81  *self = i;
82  }
83  return;
84  }
85  }
86 }
87 
88 /*--------------------------------------------------------------------------*/
89 /* vtkXMLWriterF_Delete */
90 void VTK_FORTRAN_NAME(vtkxmlwriterf_delete, VTKXMLWRITERF_DELETE)(
91  int* self
92  )
93 {
94  /* Check if the writer object exists. */
95  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
96  {
97  /* Delete this writer object. */
99 
100  /* Erase the table entry. */
101  vtkXMLWriterF_Table[*self] = 0;
102  }
103  else
104  {
105  fprintf(stderr,
106  "vtkXMLWriterF_Delete called with invalid id %d.\n",
107  *self);
108  }
109 
110  /* The writer object no longer exists. Destroy the id. */
111  *self = 0;
112 }
113 
114 /*--------------------------------------------------------------------------*/
115 /* vtkXMLWriterF_SetDataModeType */
116 void VTK_FORTRAN_NAME(vtkxmlwriterf_setdatamodetype, VTKXMLWRITERF_SETDATAMODETYPE)(
117  const int* self, const int* objType
118  )
119 {
120  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
121  {
123  }
124  else
125  {
126  fprintf(stderr,
127  "vtkXMLWriterF_SetDataModeType called with invalid id %d.\n",
128  *self);
129  }
130 }
131 
132 /*--------------------------------------------------------------------------*/
133 /* vtkXMLWriterF_SetDataObjectType */
134 void VTK_FORTRAN_NAME(vtkxmlwriterf_setdataobjecttype, VTKXMLWRITERF_SETDATAOBJECTTYPE)(
135  const int* self, const int* objType
136  )
137 {
138  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
139  {
141  }
142  else
143  {
144  fprintf(stderr,
145  "vtkXMLWriterF_SetDataObjectType called with invalid id %d.\n",
146  *self);
147  }
148 }
149 
150 /*--------------------------------------------------------------------------*/
151 /* vtkXMLWriterF_SetExtent */
152 void VTK_FORTRAN_NAME(vtkxmlwriterf_setextent, VTKXMLWRITERF_SETEXTENT)(
153  const int* self, int extent[6]
154  )
155 {
156  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
157  {
159  }
160  else
161  {
162  fprintf(stderr,
163  "vtkXMLWriterF_SetExtent called with invalid id %d.\n",
164  *self);
165  }
166 }
167 
168 /*--------------------------------------------------------------------------*/
169 /* vtkXMLWriterF_SetPoints */
170 void VTK_FORTRAN_NAME(vtkxmlwriterf_setpoints, VTKXMLWRITERF_SETPOINTS)(
171  const int* self, const int* dataType,
172  void* data, const vtkIdType* numPoints
173  )
174 {
175  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
176  {
178  data, *numPoints);
179  }
180  else
181  {
182  fprintf(stderr,
183  "vtkXMLWriterF_SetPoints called with invalid id %d.\n",
184  *self);
185  }
186 }
187 
188 /*--------------------------------------------------------------------------*/
189 /* vtkXMLWriterF_SetOrigin */
190 void VTK_FORTRAN_NAME(vtkxmlwriterf_setorigin, VTKXMLWRITERF_SETORIGIN)(
191  const int* self, double origin[3]
192  )
193 {
194  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
195  {
197  }
198  else
199  {
200  fprintf(stderr,
201  "vtkXMLWriterF_SetOrigin called with invalid id %d.\n",
202  *self);
203  }
204 }
205 
206 /*--------------------------------------------------------------------------*/
207 /* vtkXMLWriterF_SetSpacing */
208 void VTK_FORTRAN_NAME(vtkxmlwriterf_setspacing, VTKXMLWRITERF_SETSPACING)(
209  const int* self, double spacing[3]
210  )
211 {
212  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
213  {
215  }
216  else
217  {
218  fprintf(stderr,
219  "vtkXMLWriterF_SetSpacing called with invalid id %d.\n",
220  *self);
221  }
222 }
223 
224 /*--------------------------------------------------------------------------*/
225 /* vtkXMLWriterF_SetCoordinates */
226 void VTK_FORTRAN_NAME(vtkxmlwriterf_setcoordinates, VTKXMLWRITERF_SETCOORDINATES)(
227  const int* self, const int* axis, const int* dataType, void* data,
228  const vtkIdType* numCoordinates
229  )
230 {
231  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
232  {
234  *dataType, data, *numCoordinates);
235  }
236  else
237  {
238  fprintf(stderr,
239  "vtkXMLWriterF_SetCoordinates called with invalid id %d.\n",
240  *self);
241  }
242 }
243 
244 /*--------------------------------------------------------------------------*/
245 /* vtkXMLWriterF_SetCellsWithType */
246 void VTK_FORTRAN_NAME(vtkxmlwriterf_setcellswithtype, VTKXMLWRITERF_SETCELLSWITHTYPE)(
247  const int* self, const int* cellType, const vtkIdType* ncells,
248  vtkIdType* cells, const vtkIdType* cellsSize
249  )
250 {
251  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
252  {
254  *ncells, cells, *cellsSize);
255  }
256  else
257  {
258  fprintf(stderr,
259  "vtkXMLWriterF_SetCellsWithType called with invalid id %d.\n",
260  *self);
261  }
262 }
263 
264 /*--------------------------------------------------------------------------*/
265 /* vtkXMLWriterF_SetCellsWithTypes */
266 void VTK_FORTRAN_NAME(vtkxmlwriterf_setcellswithtypes, VTKXMLWRITERF_SETCELLSWITHTYPES)(
267  const int* self, int* cellTypes, const vtkIdType* ncells,
268  vtkIdType* cells, const vtkIdType* cellsSize
269  )
270 {
271  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
272  {
274  *ncells, cells, *cellsSize);
275  }
276  else
277  {
278  fprintf(stderr,
279  "vtkXMLWriterF_SetCellsWithTypes called with invalid id %d.\n",
280  *self);
281  }
282 }
283 
284 /*--------------------------------------------------------------------------*/
285 /* vtkXMLWriterF_SetPointData */
286 void VTK_FORTRAN_NAME(vtkxmlwriterf_setpointdata, VTKXMLWRITERF_SETPOINTDATA)(
287  const int* self, VTK_FORTRAN_ARG_STRING_POINTER(name),
288  const int* dataType, void* data, const vtkIdType* numTuples,
289  const int* numComponents, VTK_FORTRAN_ARG_STRING_POINTER(role)
292  )
293 {
294  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
295  {
296  /* Prepare NULL-terminated strings. */
297  const char* name_ptr = VTK_FORTRAN_REF_STRING_POINTER(name);
298  int name_length = VTK_FORTRAN_REF_STRING_LENGTH(name);
299  char* name_buffer = malloc(name_length+1);
300  const char* role_ptr = VTK_FORTRAN_REF_STRING_POINTER(role);
301  int role_length = VTK_FORTRAN_REF_STRING_LENGTH(role);
302  char* role_buffer = malloc(role_length+1);
303  if(!name_buffer || !role_buffer)
304  {
305  fprintf(stderr,
306  "vtkXMLWriterF_SetPointData failed to allocate name or role.\n");
307  if(name_buffer) { free(name_buffer); }
308  if(role_buffer) { free(role_buffer); }
309  return;
310  }
311  memcpy(name_buffer, name_ptr, name_length);
312  name_buffer[name_length] = 0;
313  memcpy(role_buffer, role_ptr, role_length);
314  role_buffer[role_length] = 0;
315 
316  /* Forward the call. */
318  *dataType, data, *numTuples, *numComponents,
319  role_buffer);
320 
321  /* Free the NULL-terminated strings. */
322  free(name_buffer);
323  free(role_buffer);
324  }
325  else
326  {
327  fprintf(stderr,
328  "vtkXMLWriterF_SetPointData called with invalid id %d.\n",
329  *self);
330  }
331 }
332 
333 /*--------------------------------------------------------------------------*/
334 /* vtkXMLWriterF_SetCellData */
335 void VTK_FORTRAN_NAME(vtkxmlwriterf_setcelldata, VTKXMLWRITERF_SETCELLDATA)(
336  const int* self, VTK_FORTRAN_ARG_STRING_POINTER(name),
337  const int* dataType, void* data, const vtkIdType* numTuples,
338  const int* numComponents, VTK_FORTRAN_ARG_STRING_POINTER(role)
341  )
342 {
343  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
344  {
345  /* Prepare NULL-terminated strings. */
346  const char* name_ptr = VTK_FORTRAN_REF_STRING_POINTER(name);
347  int name_length = VTK_FORTRAN_REF_STRING_LENGTH(name);
348  char* name_buffer = malloc(name_length+1);
349  const char* role_ptr = VTK_FORTRAN_REF_STRING_POINTER(role);
350  int role_length = VTK_FORTRAN_REF_STRING_LENGTH(role);
351  char* role_buffer = malloc(role_length+1);
352  if(!name_buffer || !role_buffer)
353  {
354  fprintf(stderr,
355  "vtkXMLWriterF_SetCellData failed to allocate name or role.\n");
356  if(name_buffer) { free(name_buffer); }
357  if(role_buffer) { free(role_buffer); }
358  return;
359  }
360  memcpy(name_buffer, name_ptr, name_length);
361  name_buffer[name_length] = 0;
362  memcpy(role_buffer, role_ptr, role_length);
363  role_buffer[role_length] = 0;
364 
365  /* Forward the call. */
367  *dataType, data, *numTuples, *numComponents,
368  role_buffer);
369 
370  /* Free the NULL-terminated strings. */
371  free(name_buffer);
372  free(role_buffer);
373  }
374  else
375  {
376  fprintf(stderr,
377  "vtkXMLWriterF_SetCellData called with invalid id %d.\n",
378  *self);
379  }
380 }
381 
382 /*--------------------------------------------------------------------------*/
383 /* vtkXMLWriterF_SetFileName */
384 void VTK_FORTRAN_NAME(vtkxmlwriterf_setfilename, VTKXMLWRITERF_SETFILENAME)(
385  const int* self, VTK_FORTRAN_ARG_STRING_POINTER(name)
387  )
388 {
389  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
390  {
391  /* Prepare NULL-terminated string. */
392  const char* name_ptr = VTK_FORTRAN_REF_STRING_POINTER(name);
393  int name_length = VTK_FORTRAN_REF_STRING_LENGTH(name);
394  char* name_buffer = malloc(name_length+1);
395  if(!name_buffer)
396  {
397  fprintf(stderr,
398  "vtkXMLWriterF_SetFileName failed to allocate name.\n");
399  return;
400  }
401  memcpy(name_buffer, name_ptr, name_length);
402  name_buffer[name_length] = 0;
403 
404  /* Forward the call. */
405  vtkXMLWriterC_SetFileName(vtkXMLWriterF_Table[*self], name_buffer);
406 
407  /* Free the NULL-terminated string. */
408  free(name_buffer);
409  }
410  else
411  {
412  fprintf(stderr,
413  "vtkXMLWriterF_SetFileName called with invalid id %d.\n",
414  *self);
415  }
416 }
417 
418 /*--------------------------------------------------------------------------*/
419 /* vtkXMLWriterF_Write */
420 void VTK_FORTRAN_NAME(vtkxmlwriterf_write, VTKXMLWRITERF_WRITE)(
421  const int* self, int* success
422  )
423 {
424  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
425  {
426  *success = vtkXMLWriterC_Write(vtkXMLWriterF_Table[*self]);
427  }
428  else
429  {
430  fprintf(stderr,
431  "vtkXMLWriterF_Write called with invalid id %d.\n",
432  *self);
433  }
434 }
435 
436 /*--------------------------------------------------------------------------*/
437 /* vtkXMLWriterF_SetNumberOfTimeSteps */
438 void VTK_FORTRAN_NAME(vtkxmlwriterf_setnumberoftimesteps, VTKXMLWRITERF_SETNUMBEROFTIMESTEPS)(
439  const int* self, const int* numTimeSteps
440  )
441 {
442  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
443  {
445  *numTimeSteps);
446  }
447  else
448  {
449  fprintf(stderr,
450  "vtkXMLWriterF_SetNumberOfTimeSteps called with invalid id %d.\n",
451  *self);
452  }
453 }
454 
455 /*--------------------------------------------------------------------------*/
456 /* vtkXMLWriterF_Start */
457 void VTK_FORTRAN_NAME(vtkxmlwriterf_start, VTKXMLWRITERF_START)(
458  const int* self
459  )
460 {
461  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
462  {
464  }
465  else
466  {
467  fprintf(stderr,
468  "vtkXMLWriterF_Start called with invalid id %d.\n",
469  *self);
470  }
471 }
472 
473 /*--------------------------------------------------------------------------*/
474 /* vtkXMLWriterF_WriteNextTimeStep */
475 void VTK_FORTRAN_NAME(vtkxmlwriterf_writenexttimestep, VTKXMLWRITERF_WRITENEXTTIMESTEP)(
476  const int* self, const double* timeValue
477  )
478 {
479  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
480  {
482  }
483  else
484  {
485  fprintf(stderr,
486  "vtkXMLWriterF_WriteNextTimeStep called with invalid id %d.\n",
487  *self);
488  }
489 }
490 
491 /*--------------------------------------------------------------------------*/
492 /* vtkXMLWriterF_Stop */
493 void VTK_FORTRAN_NAME(vtkxmlwriterf_stop, VTKXMLWRITERF_STOP)(
494  const int* self
495  )
496 {
497  if(*self > 0 && *self <= VTK_XMLWRITERF_MAX && vtkXMLWriterF_Table[*self])
498  {
500  }
501  else
502  {
503  fprintf(stderr,
504  "vtkXMLWriterF_Stop called with invalid id %d.\n",
505  *self);
506  }
507 }
508 // VTK-HeaderTest-Exclude: vtkXMLWriterF.h
struct vtkXMLWriterC_s vtkXMLWriterC
vtkXMLWriterC is an opaque structure holding the state of an individual writer object.
Definition: vtkXMLWriterC.h:30
#define VTK_FORTRAN_REF_STRING_LENGTH(name)
Definition: vtkXMLWriterF.h:59
VTKIOXML_EXPORT void vtkXMLWriterC_SetCellsWithTypes(vtkXMLWriterC *self, int *cellTypes, vtkIdType ncells, vtkIdType *cells, vtkIdType cellsSize)
Set a cell array on the data object to be written.
VTKIOXML_EXPORT void vtkXMLWriterC_SetPointData(vtkXMLWriterC *self, const char *name, int dataType, void *data, vtkIdType numTuples, int numComponents, const char *role)
Set a point or cell data array by name.
VTKIOXML_EXPORT void vtkXMLWriterC_SetCellsWithType(vtkXMLWriterC *self, int cellType, vtkIdType ncells, vtkIdType *cells, vtkIdType cellsSize)
Set a cell array on the data object to be written.
#define VTK_FORTRAN_ARG_STRING_POINTER(name)
Definition: vtkXMLWriterF.h:50
VTKIOXML_EXPORT void vtkXMLWriterC_Delete(vtkXMLWriterC *self)
Delete the writer object.
#define VTK_XMLWRITERF_MAX
Definition: vtkXMLWriterF.h:37
VTKIOXML_EXPORT void vtkXMLWriterC_SetPoints(vtkXMLWriterC *self, int dataType, void *data, vtkIdType numPoints)
Set the points of a point data set.
#define VTK_FORTRAN_NAME(name, NAME)
Definition: vtkXMLWriterF.h:47
int vtkIdType
Definition: vtkType.h:275
VTKIOXML_EXPORT void vtkXMLWriterC_Start(vtkXMLWriterC *self)
Start writing a time-series to the output file.
VTKIOXML_EXPORT void vtkXMLWriterC_SetNumberOfTimeSteps(vtkXMLWriterC *self, int numTimeSteps)
Set the number of time steps that will be written between upcoming Start and Stop calls...
VTKIOXML_EXPORT void vtkXMLWriterC_SetOrigin(vtkXMLWriterC *self, double origin[3])
Set the origin of an image data set.
VTKIOXML_EXPORT void vtkXMLWriterC_SetDataObjectType(vtkXMLWriterC *self, int objType)
Set the VTK data object type that will be written.
VTKIOXML_EXPORT void vtkXMLWriterC_SetCoordinates(vtkXMLWriterC *self, int axis, int dataType, void *data, vtkIdType numCoordinates)
Set the coordinates along one axis of a rectilinear grid data set.
#define VTK_FORTRAN_ARG_STRING_LENGTH(name)
Definition: vtkXMLWriterF.h:53
VTKIOXML_EXPORT void vtkXMLWriterC_SetDataModeType(vtkXMLWriterC *self, int datamodetype)
Set the VTK writer data mode to either:
VTKIOXML_EXPORT int vtkXMLWriterC_Write(vtkXMLWriterC *self)
Write the data to a file immediately.
VTKIOXML_EXPORT void vtkXMLWriterC_SetSpacing(vtkXMLWriterC *self, double spacing[3])
Set the spacing of an image data set.
VTKIOXML_EXPORT void vtkXMLWriterC_SetFileName(vtkXMLWriterC *self, const char *fileName)
Set the name of the file into which the data are to be written.
#define VTK_FORTRAN_REF_STRING_POINTER(name)
Definition: vtkXMLWriterF.h:56
static vtkXMLWriterC * vtkXMLWriterF_Table[VTK_XMLWRITERF_MAX+1]
Definition: vtkXMLWriterF.h:38
VTKIOXML_EXPORT vtkXMLWriterC * vtkXMLWriterC_New()
Create a new instance of vtkXMLWriterC.
CellTypeInDataSet cellType(vtkDataSet *input)
VTKIOXML_EXPORT void vtkXMLWriterC_WriteNextTimeStep(vtkXMLWriterC *self, double timeValue)
Write one time step of a time-series to the output file.
VTKIOXML_EXPORT void vtkXMLWriterC_SetExtent(vtkXMLWriterC *self, int extent[6])
Set the extent of a structured data set.
VTKIOXML_EXPORT void vtkXMLWriterC_SetCellData(vtkXMLWriterC *self, const char *name, int dataType, void *data, vtkIdType numTuples, int numComponents, const char *role)
VTKIOXML_EXPORT void vtkXMLWriterC_Stop(vtkXMLWriterC *self)
Stop writing a time-series to the output file.