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