hydrocode_Radial_Lag 0.3
This is an implementation of fully explict forward Euler scheme for multi-D radially symmetric compressible flows on Lagrangian coordinate
file_out_hdf5.c 文件参考

This is a set of functions which control the readout of 1-D and 2-D data in HDF5 file format. 更多...

#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include "../include/var_struc.h"
#include "../include/file_io.h"
#include "hdf5.h"
file_out_hdf5.c 的引用(Include)关系图:

浏览源代码.

宏定义

#define PRINT_NC(v, v_array)
 Print out fluid variable 'v' with data array 'v_array'. 更多...
 

函数

void file_1D_write_HDF5 (const int m, const int N, const struct cell_var_stru CV, double *X[], const double *cpu_time, const char *problem, double time_plot[])
 This function write the 1-D solution into HDF5 output '.h5' files. 更多...
 
void file_2D_write_HDF5 (const int n_x, const int n_y, const int N, const struct cell_var_stru CV[], double **X, double **Y, const double *cpu_time, const char *problem, double time_plot[])
 This function write the 2-D solution into HDF5 output '.h5' files. 更多...
 

详细描述

This is a set of functions which control the readout of 1-D and 2-D data in HDF5 file format.

注意
Library Dependency: HDF5®

在文件 file_out_hdf5.c 中定义.

宏定义说明

◆ PRINT_NC

#define PRINT_NC (   v,
  v_array 
)
值:
do { \
dataset_id = H5Dcreate(group_id, #v, H5T_NATIVE_FLOAT, dataspace_id, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT); \
status = H5Dwrite(dataset_id, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, v_array); \
status = H5Dclose(dataset_id); \
} while (0)

Print out fluid variable 'v' with data array 'v_array'.

在文件 file_out_hdf5.c37 行定义.

函数说明

◆ file_1D_write_HDF5()

void file_1D_write_HDF5 ( const int  m,
const int  N,
const struct cell_var_stru  CV,
double *  X[],
const double *  cpu_time,
const char *  problem,
double  time_plot[] 
)

This function write the 1-D solution into HDF5 output '.h5' files.

参数
[in]mThe number of spatial points in the output data.
[in]NThe number of time steps in the output data.
[in]CVStructure of grid variable data in computational grid cells.
[in]X[]Array of the coordinate data.
[in]cpu_timeArray of the CPU time recording.
[in]problemName of the numerical results for the test problem.
[in]time_plotArray of the plotting time recording.

在文件 file_out_hdf5.c54 行定义.

函数调用图:
这是这个函数的调用关系图:

◆ file_2D_write_HDF5()

void file_2D_write_HDF5 ( const int  n_x,
const int  n_y,
const int  N,
const struct cell_var_stru  CV[],
double **  X,
double **  Y,
const double *  cpu_time,
const char *  problem,
double  time_plot[] 
)

This function write the 2-D solution into HDF5 output '.h5' files.

参数
[in]n_xThe number of x-spatial points in the output data.
[in]n_yThe number of y-spatial points in the output data.
[in]NThe number of time steps in the output data.
[in]CVStructure of variable data in computational grid cells.
[in]XArray of the x-coordinate data.
[in]YArray of the y-coordinate data.
[in]cpu_timeArray of the CPU time recording.
[in]problemName of the numerical results for the test problem.
[in]time_plotArray of the plotting time recording.

在文件 file_out_hdf5.c162 行定义.

函数调用图: