SDDSlib
Loading...
Searching...
No Matches
SDDSmpi_input.c File Reference

SDDS MPI Input Initialization and Data Broadcasting Functions. More...

#include "SDDS.h"
#include "SDDS_internal.h"
#include "scan.h"

Go to the source code of this file.

Classes

struct  ELEMENT_DEF
 Structure defining an element with various attributes. More...
 
struct  ASSOCIATE_DEF
 Structure defining an associate with various attributes. More...
 
struct  OTHER_DEF
 Structure defining additional layout information. More...
 
struct  STRING_DEF
 Structure defining a string with fixed maximum length. More...
 

Macros

#define STRING_COL_LENGTH   40
 Defines the maximum length of strings in arrays or columns.
 

Functions

int32_t SDDS_MPI_ReadPage (SDDS_DATASET *SDDS_dataset)
 Reads a page from an SDDS dataset using MPI.
 
int32_t SDDS_MPI_BroadcastLayout (SDDS_DATASET *SDDS_dataset)
 Broadcasts the layout of an SDDS dataset to all MPI processes.
 
int32_t SDDS_MPI_InitializeInput (SDDS_DATASET *SDDS_dataset, char *filename)
 Initializes an SDDS dataset for input using MPI.
 
int32_t SDDS_MPI_InitializeInputFromSearchPath (SDDS_DATASET *SDDS_dataset, char *file)
 Initializes an SDDS dataset for input by searching the provided search path using MPI.
 
int32_t SDDS_Master_InitializeInputFromSearchPath (SDDS_DATASET *SDDS_dataset, MPI_DATASET *MPI_dataset, char *file)
 Initializes an SDDS dataset for input from a search path on the master MPI process.
 
int32_t SDDS_Master_InitializeInput (SDDS_DATASET *SDDS_dataset, MPI_DATASET *MPI_dataset, char *file)
 Initializes an SDDS dataset for input on the master MPI process.
 
int32_t SDDS_Master_ReadPage (SDDS_DATASET *SDDS_dataset)
 Reads a single page of data from an SDDS dataset on the master MPI process and broadcasts it.
 

Detailed Description

SDDS MPI Input Initialization and Data Broadcasting Functions.

This file implements the core functionalities for initializing and reading Self Describing Data Sets (SDDS) in a parallel computing environment using MPI. It includes definitions of essential data structures and functions that handle various data formats, manage layout information, and facilitate the broadcasting of data across multiple MPI processes.

License
This file is distributed under the terms of the Software License Agreement found in the file LICENSE included with this distribution.
Authors
H. Shang R. Soliday

Definition in file SDDSmpi_input.c.

Macro Definition Documentation

◆ STRING_COL_LENGTH

#define STRING_COL_LENGTH   40

Defines the maximum length of strings in arrays or columns.

The string length in arrays or columns should be less than 40 characters.

Definition at line 95 of file SDDSmpi_input.c.

Function Documentation

◆ SDDS_Master_InitializeInput()

int32_t SDDS_Master_InitializeInput ( SDDS_DATASET * SDDS_dataset,
MPI_DATASET * MPI_dataset,
char * file )

Initializes an SDDS dataset for input on the master MPI process.

This function initializes the SDDS_DATASET structure for reading from a specified file using MPI. It is intended to be called by the master (root) process, which performs the dataset initialization and then broadcasts the layout to all other MPI processes. Non-master processes prepare their datasets based on the received layout information.

Parameters
[in,out]SDDS_datasetPointer to the SDDS_DATASET structure to be initialized.
[in]MPI_datasetPointer to the MPI_DATASET structure containing MPI-related information such as the communicator and process ID.
[in]fileThe name of the file to be opened and initialized for reading.
Returns
  • 1 on successful initialization and layout broadcast.
  • 0 if initialization fails or layout broadcasting fails.

The function performs the following steps:

  • On the master process (MPI rank 0), it initializes the dataset using SDDS_InitializeInput.
  • Non-master processes zero out their SDDS_DATASET structures to prepare for receiving the layout.
  • Sets the parallel_io flag to 0 and assigns the MPI_dataset.
  • Broadcasts the layout information to all MPI processes using SDDS_MPI_BroadcastLayout.
Note
  • Only the master process performs dataset initialization.
  • Non-master processes rely entirely on the broadcasted layout information.
See also
SDDS_MPI_BroadcastLayout, SDDS_InitializeInput, SDDS_ZeroMemory

Definition at line 877 of file SDDSmpi_input.c.

877 {
878 if (MPI_dataset->myid == 0) {
879 if (!SDDS_InitializeInput(SDDS_dataset, file))
880 return 0;
881 } else {
882 if (!SDDS_ZeroMemory((void *)SDDS_dataset, sizeof(SDDS_DATASET)))
883 SDDS_Bomb("Unable to zero memory for SDDS dataset(SDDS_MPI_Setup)");
884 }
885 SDDS_dataset->parallel_io = 0;
886 SDDS_dataset->MPI_dataset = MPI_dataset;
887 if (!SDDS_MPI_BroadcastLayout(SDDS_dataset))
888 return 0;
889
890 return 1;
891}
int32_t SDDS_InitializeInput(SDDS_DATASET *SDDS_dataset, char *filename)
Definition SDDS_input.c:49
int32_t SDDS_ZeroMemory(void *mem, int64_t n_bytes)
Sets a block of memory to zero.
void SDDS_Bomb(char *message)
Terminates the program after printing an error message and recorded errors.
Definition SDDS_utils.c:342
int32_t SDDS_MPI_BroadcastLayout(SDDS_DATASET *SDDS_dataset)
Broadcasts the layout of an SDDS dataset to all MPI processes.

◆ SDDS_Master_InitializeInputFromSearchPath()

int32_t SDDS_Master_InitializeInputFromSearchPath ( SDDS_DATASET * SDDS_dataset,
MPI_DATASET * MPI_dataset,
char * file )

Initializes an SDDS dataset for input from a search path on the master MPI process.

This function is designed to be called by the master (root) process in an MPI environment. It searches for the specified file within the search paths, initializes the SDDS_DATASET structure, and broadcasts the layout information to all other MPI processes. Non-master processes receive the broadcasted layout and prepare their datasets accordingly.

Parameters
[in,out]SDDS_datasetPointer to the SDDS_DATASET structure to be initialized.
[in]MPI_datasetPointer to the MPI_DATASET structure containing MPI-related information such as the communicator and process ID.
[in]fileThe name of the file to search for and initialize. The function searches for this file in the configured search paths.
Returns
  • 1 on successful initialization and layout broadcast.
  • 0 if the file is not found, initialization fails, or layout broadcasting fails.

The function performs the following steps:

  • On the master process (MPI rank 0), it searches for the specified file in the search paths.
  • If the file is found, it initializes the dataset using SDDS_InitializeInput.
  • Non-master processes zero out their SDDS_DATASET structures to prepare for receiving the layout.
  • Sets the parallel_io flag to 0 and assigns the MPI_dataset.
  • Broadcasts the layout information to all MPI processes using SDDS_MPI_BroadcastLayout.
Note
  • Only the master process performs file searching and dataset initialization.
  • Non-master processes rely entirely on the broadcasted layout information.
See also
SDDS_MPI_InitializeInputFromSearchPath, SDDS_MPI_BroadcastLayout, SDDS_InitializeInput, SDDS_ZeroMemory

Definition at line 816 of file SDDSmpi_input.c.

816 {
817 char *filename;
818
819 if (MPI_dataset->myid == 0) {
820 if (!(filename = findFileInSearchPath(file))) {
821 SDDS_SetError("file does not exist in searchpath (InitializeInputFromSearchPath)");
822 return 0;
823 }
824 if (!SDDS_InitializeInput(SDDS_dataset, filename)) {
825 free(filename);
826 return 0;
827 }
828 free(filename);
829 } else {
830 if (!SDDS_ZeroMemory((void *)SDDS_dataset, sizeof(SDDS_DATASET)))
831 SDDS_Bomb("Unable to zero memory for SDDS dataset(SDDS_MPI_Setup)");
832 }
833 SDDS_dataset->parallel_io = 0;
834 SDDS_dataset->MPI_dataset = MPI_dataset;
835 if (!SDDS_MPI_BroadcastLayout(SDDS_dataset))
836 return 0;
837
838 return 1;
839}
void SDDS_SetError(char *error_text)
Records an error message in the SDDS error stack.
Definition SDDS_utils.c:379

◆ SDDS_Master_ReadPage()

int32_t SDDS_Master_ReadPage ( SDDS_DATASET * SDDS_dataset)

Reads a single page of data from an SDDS dataset on the master MPI process and broadcasts it.

This function is intended to be called by the master (root) MPI process. It reads a single page of data from the SDDS input file, counts the number of rows of interest, and then broadcasts the data to all other MPI processes. Non-master processes allocate memory for the received data and populate their datasets accordingly.

Parameters
[in,out]SDDS_datasetPointer to the SDDS_DATASET structure from which the page will be read and into which the data will be stored.
Returns
  • The number of rows read on success.
  • A negative value if an error occurs during reading or broadcasting.

The function performs the following steps:

  • On the master process (MPI rank 0), it reads a page of data using SDDS_ReadPage.
  • Counts the number of rows of interest using SDDS_CountRowsOfInterest.
  • Broadcasts the number of rows and retrieval status to all MPI processes.
  • Non-master processes allocate memory for the received data by starting a new page.
  • Broadcasts each parameter's value to all processes based on its type.
  • Broadcasts each array's data to all processes based on its type.
  • Broadcasts the data for each column to all processes based on its type.
  • Allocates and populates string data where necessary.
Note
  • The function assumes that only one page of data is present in the input file.
  • String data requires special handling to allocate appropriate memory on non-master processes.
See also
SDDS_ReadPage, SDDS_CountRowsOfInterest, SDDS_StartPage, SDDS_DefineColumn, SDDS_DefineParameter, SDDS_DefineArray, SDDS_DefineAssociate

Definition at line 928 of file SDDSmpi_input.c.

928 {
929 MPI_DATASET *MPI_dataset = SDDS_dataset->MPI_dataset;
930 int32_t rows = 0, i, j, len, retrival = 0;
931 STRING_DEF *str_val = NULL;
932
933 /*master read file */
934 if (MPI_dataset->myid == 0) {
935 if ((retrival = SDDS_ReadPage(SDDS_dataset)) <= 0) {
936 SDDS_SetError("SDDS_MPI_ReadParameterFile2: Error in reading input file");
937 return (retrival);
938 }
939 rows = SDDS_CountRowsOfInterest(SDDS_dataset);
940 }
941 MPI_Bcast(&rows, 1, MPI_INT, 0, MPI_dataset->comm);
942 MPI_Bcast(&retrival, 1, MPI_INT, 0, MPI_dataset->comm);
943 if (MPI_dataset->myid != 0) {
944 /*allocate memory for other processors */
945 if (!SDDS_StartPage(SDDS_dataset, rows))
946 return 0;
947 }
948 /*broadcast the parameters to other processors */
949 for (i = 0; i < SDDS_dataset->layout.n_parameters; i++) {
950 switch (SDDS_dataset->layout.parameter_definition[i].type) {
951 case SDDS_LONGDOUBLE:
952 MPI_Bcast(SDDS_dataset->parameter[i], 1, MPI_LONG_DOUBLE, 0, MPI_dataset->comm);
953 break;
954 case SDDS_DOUBLE:
955 MPI_Bcast(SDDS_dataset->parameter[i], 1, MPI_DOUBLE, 0, MPI_dataset->comm);
956 break;
957 case SDDS_FLOAT:
958 MPI_Bcast(SDDS_dataset->parameter[i], 1, MPI_FLOAT, 0, MPI_dataset->comm);
959 break;
960 case SDDS_LONG64:
961 MPI_Bcast(SDDS_dataset->parameter[i], 1, MPI_INT64_T, 0, MPI_dataset->comm);
962 break;
963 case SDDS_ULONG64:
964 MPI_Bcast(SDDS_dataset->parameter[i], 1, MPI_UINT64_T, 0, MPI_dataset->comm);
965 break;
966 case SDDS_LONG:
967 MPI_Bcast(SDDS_dataset->parameter[i], 1, MPI_INT, 0, MPI_dataset->comm);
968 break;
969 case SDDS_ULONG:
970 MPI_Bcast(SDDS_dataset->parameter[i], 1, MPI_UNSIGNED, 0, MPI_dataset->comm);
971 break;
972 case SDDS_SHORT:
973 MPI_Bcast(SDDS_dataset->parameter[i], 1, MPI_SHORT, 0, MPI_dataset->comm);
974 break;
975 case SDDS_USHORT:
976 MPI_Bcast(SDDS_dataset->parameter[i], 1, MPI_UNSIGNED_SHORT, 0, MPI_dataset->comm);
977 break;
978 case SDDS_STRING:
979 if (MPI_dataset->myid == 0)
980 len = strlen(*(char **)SDDS_dataset->parameter[i]);
981 MPI_Bcast(&len, 1, MPI_INT, 0, MPI_dataset->comm);
982 if (MPI_dataset->myid != 0)
983 SDDS_dataset->parameter[i] = (char *)malloc(sizeof(char) * len);
984 MPI_Bcast(SDDS_dataset->parameter[i], len, MPI_BYTE, 0, MPI_dataset->comm);
985 break;
986 case SDDS_CHARACTER:
987 MPI_Bcast(SDDS_dataset->parameter[i], 1, MPI_CHAR, 0, MPI_dataset->comm);
988 break;
989 }
990 }
991 /* broadcast arrays to other processros */
992 for (i = 0; i < SDDS_dataset->layout.n_arrays; i++) {
993 MPI_Bcast(&(SDDS_dataset->layout.array_definition[i].dimensions), 1, MPI_INT, 0, MPI_dataset->comm);
994 switch (SDDS_dataset->layout.array_definition[i].type) {
995 case SDDS_LONGDOUBLE:
996 MPI_Bcast(SDDS_dataset->array[i].data, SDDS_dataset->layout.array_definition[i].dimensions, MPI_LONG_DOUBLE, 0, MPI_dataset->comm);
997 break;
998 case SDDS_DOUBLE:
999 MPI_Bcast(SDDS_dataset->array[i].data, SDDS_dataset->layout.array_definition[i].dimensions, MPI_DOUBLE, 0, MPI_dataset->comm);
1000 break;
1001 case SDDS_FLOAT:
1002 MPI_Bcast(SDDS_dataset->array[i].data, SDDS_dataset->layout.array_definition[i].dimensions, MPI_FLOAT, 0, MPI_dataset->comm);
1003 break;
1004 case SDDS_LONG64:
1005 MPI_Bcast(SDDS_dataset->array[i].data, SDDS_dataset->layout.array_definition[i].dimensions, MPI_INT64_T, 0, MPI_dataset->comm);
1006 break;
1007 case SDDS_ULONG64:
1008 MPI_Bcast(SDDS_dataset->array[i].data, SDDS_dataset->layout.array_definition[i].dimensions, MPI_UINT64_T, 0, MPI_dataset->comm);
1009 break;
1010 case SDDS_LONG:
1011 MPI_Bcast(SDDS_dataset->array[i].data, SDDS_dataset->layout.array_definition[i].dimensions, MPI_INT, 0, MPI_dataset->comm);
1012 break;
1013 case SDDS_ULONG:
1014 MPI_Bcast(SDDS_dataset->array[i].data, SDDS_dataset->layout.array_definition[i].dimensions, MPI_UNSIGNED, 0, MPI_dataset->comm);
1015 break;
1016 case SDDS_SHORT:
1017 MPI_Bcast(SDDS_dataset->array[i].data, SDDS_dataset->layout.array_definition[i].dimensions, MPI_SHORT, 0, MPI_dataset->comm);
1018 break;
1019 case SDDS_USHORT:
1020 MPI_Bcast(SDDS_dataset->array[i].data, SDDS_dataset->layout.array_definition[i].dimensions, MPI_UNSIGNED_SHORT, 0, MPI_dataset->comm);
1021 break;
1022 case SDDS_STRING:
1023 str_val = malloc(sizeof(*str_val) * SDDS_dataset->layout.array_definition[i].dimensions);
1024 if (MPI_dataset->myid == 0) {
1025 for (j = 0; j < SDDS_dataset->layout.array_definition[i].dimensions; j++)
1026 sprintf(str_val[j].str_value, "%s", ((char **)(SDDS_dataset->array[i].data))[j]);
1027 }
1028 MPI_Bcast(str_val, SDDS_dataset->layout.array_definition[i].dimensions * 256, MPI_BYTE, 0, MPI_dataset->comm);
1029 if (MPI_dataset->myid) {
1030 for (j = 0; j < SDDS_dataset->layout.array_definition[i].dimensions; j++)
1031 SDDS_CopyString(&(((char **)(SDDS_dataset->array[i].data))[j]), str_val[j].str_value);
1032 }
1033 free(str_val);
1034 str_val = NULL;
1035 break;
1036 case SDDS_CHARACTER:
1037 MPI_Bcast(SDDS_dataset->array[i].data, SDDS_dataset->layout.array_definition[i].dimensions, MPI_CHAR, 0, MPI_dataset->comm);
1038 break;
1039 }
1040 }
1041
1042 /* broadcast columns to other processors */
1043 SDDS_dataset->n_rows = SDDS_dataset->n_rows_allocated = rows;
1044 MPI_Bcast(SDDS_dataset->row_flag, rows, MPI_INT, 0, MPI_dataset->comm);
1045 for (i = 0; i < SDDS_dataset->layout.n_columns; i++) {
1046 switch (SDDS_dataset->layout.column_definition[i].type) {
1047 case SDDS_LONGDOUBLE:
1048 MPI_Bcast(SDDS_dataset->data[i], rows, MPI_LONG_DOUBLE, 0, MPI_dataset->comm);
1049 break;
1050 case SDDS_DOUBLE:
1051 MPI_Bcast(SDDS_dataset->data[i], rows, MPI_DOUBLE, 0, MPI_dataset->comm);
1052 break;
1053 case SDDS_FLOAT:
1054 MPI_Bcast(SDDS_dataset->data[i], rows, MPI_FLOAT, 0, MPI_dataset->comm);
1055 break;
1056 case SDDS_LONG64:
1057 MPI_Bcast(SDDS_dataset->data[i], rows, MPI_INT64_T, 0, MPI_dataset->comm);
1058 break;
1059 case SDDS_ULONG64:
1060 MPI_Bcast(SDDS_dataset->data[i], rows, MPI_UINT64_T, 0, MPI_dataset->comm);
1061 break;
1062 case SDDS_LONG:
1063 MPI_Bcast(SDDS_dataset->data[i], rows, MPI_INT, 0, MPI_dataset->comm);
1064 break;
1065 case SDDS_ULONG:
1066 MPI_Bcast(SDDS_dataset->data[i], rows, MPI_UNSIGNED, 0, MPI_dataset->comm);
1067 break;
1068 case SDDS_SHORT:
1069 MPI_Bcast(SDDS_dataset->data[i], rows, MPI_SHORT, 0, MPI_dataset->comm);
1070 break;
1071 case SDDS_USHORT:
1072 MPI_Bcast(SDDS_dataset->data[i], rows, MPI_UNSIGNED_SHORT, 0, MPI_dataset->comm);
1073 break;
1074 case SDDS_STRING:
1075 if (!str_val)
1076 str_val = malloc(sizeof(*str_val) * rows);
1077 if (MPI_dataset->myid == 0) {
1078 for (j = 0; j < rows; j++)
1079 sprintf(str_val[j].str_value, "%s", ((char **)SDDS_dataset->data[i])[j]);
1080 }
1081 MPI_Bcast(str_val, rows * STRING_COL_LENGTH, MPI_BYTE, 0, MPI_dataset->comm);
1082 if (MPI_dataset->myid) {
1083 for (j = 0; j < rows; j++)
1084 SDDS_CopyString(&(((char **)SDDS_dataset->data[i])[j]), str_val[j].str_value);
1085 }
1086 break;
1087 case SDDS_CHARACTER:
1088 MPI_Bcast(SDDS_dataset->data[i], rows, MPI_CHAR, 0, MPI_dataset->comm);
1089 break;
1090 }
1091 }
1092 if (str_val)
1093 free(str_val);
1094 return retrival;
1095}
int32_t SDDS_StartPage(SDDS_DATASET *SDDS_dataset, int64_t expected_n_rows)
int64_t SDDS_CountRowsOfInterest(SDDS_DATASET *SDDS_dataset)
Counts the number of rows marked as "of interest" in the current data table.
int32_t SDDS_ReadPage(SDDS_DATASET *SDDS_dataset)
int32_t SDDS_CopyString(char **target, const char *source)
Copies a source string to a target string with memory allocation.
Definition SDDS_utils.c:856
#define STRING_COL_LENGTH
Defines the maximum length of strings in arrays or columns.
#define SDDS_ULONG
Identifier for the unsigned 32-bit integer data type.
Definition SDDStypes.h:67
#define SDDS_FLOAT
Identifier for the float data type.
Definition SDDStypes.h:43
#define SDDS_STRING
Identifier for the string data type.
Definition SDDStypes.h:85
#define SDDS_ULONG64
Identifier for the unsigned 64-bit integer data type.
Definition SDDStypes.h:55
#define SDDS_LONG
Identifier for the signed 32-bit integer data type.
Definition SDDStypes.h:61
#define SDDS_SHORT
Identifier for the signed short integer data type.
Definition SDDStypes.h:73
#define SDDS_CHARACTER
Identifier for the character data type.
Definition SDDStypes.h:91
#define SDDS_USHORT
Identifier for the unsigned short integer data type.
Definition SDDStypes.h:79
#define SDDS_DOUBLE
Identifier for the double data type.
Definition SDDStypes.h:37
#define SDDS_LONGDOUBLE
Identifier for the long double data type.
Definition SDDStypes.h:31
#define SDDS_LONG64
Identifier for the signed 64-bit integer data type.
Definition SDDStypes.h:49
Structure defining a string with fixed maximum length.

◆ SDDS_MPI_BroadcastLayout()

int32_t SDDS_MPI_BroadcastLayout ( SDDS_DATASET * SDDS_dataset)

Broadcasts the layout of an SDDS dataset to all MPI processes.

This function creates MPI data types for the layout structures, broadcasts the layout information from the root process to all other processes, and defines columns, parameters, arrays, and associates on non-root processes based on the received layout information. It ensures that all MPI processes have a consistent view of the dataset layout.

Parameters
SDDS_datasetPointer to the SDDS_DATASET structure representing the dataset.
Returns
  • 1 on successful broadcast and layout definition.
  • 0 if an error occurs during memory allocation, MPI operations, or layout definition.

The function performs the following steps:

  • Checks and commits MPI data types for ELEMENT_DEF and OTHER_DEF structures.
  • On the root process (MPI rank 0), fills the OTHER_DEF structure with layout information.
  • Broadcasts the OTHER_DEF structure to all processes.
  • Allocates memory for columns, parameters, arrays, and associates based on the received counts.
  • On the root process, populates the ELEMENT_DEF and ASSOCIATE_DEF arrays with layout details.
  • Broadcasts the ELEMENT_DEF and ASSOCIATE_DEF arrays to all processes.
  • On non-root processes, defines columns, parameters, arrays, and associates using the received information.
  • Saves the layout on non-root processes.
See also
SDDS_DefineColumn, SDDS_DefineParameter, SDDS_DefineArray, SDDS_DefineAssociate, SDDS_SaveLayout

Definition at line 171 of file SDDSmpi_input.c.

171 {
172 SDDS_LAYOUT *layout = &(SDDS_dataset->layout);
173 MPI_DATASET *MPI_dataset = SDDS_dataset->MPI_dataset;
174 char *symbol, *units, *description, *format_string, *fixed_value, *filename, *path, *contents;
175 int i;
176 ELEMENT_DEF *column = NULL, *parameter = NULL, *array = NULL;
177 ASSOCIATE_DEF *associate = NULL;
178 OTHER_DEF other;
179 MPI_Datatype elementType, otherType, oldtypes[4], associateType;
180 int blockcounts[4];
181 /* MPI_Aint type used to be consistent with syntax of */
182 /* MPI_Type_extent routine */
183 MPI_Aint offsets[4], int_ext, short_ext, uint_ext;
184 MPI_Aint int_lb, short_lb, uint_lb;
185
186 MPI_Type_get_extent(MPI_INT, &int_lb, &int_ext);
187 MPI_Type_get_extent(MPI_UNSIGNED, &uint_lb, &uint_ext);
188 MPI_Type_get_extent(MPI_SHORT, &short_lb, &short_ext);
189
190 layout = &(SDDS_dataset->layout);
191 /*commit element type */
192 offsets[0] = 0;
193 oldtypes[0] = MPI_INT;
194 blockcounts[0] = 13;
195 offsets[1] = 13 * int_ext;
196 oldtypes[1] = MPI_CHAR;
197 blockcounts[1] = 256 + 256 + 256 + 1024 + 256 + 1024 + 256;
198 /* Now define structured type and commit it */
199 MPI_Type_create_struct(2, blockcounts, offsets, oldtypes, &elementType);
200 MPI_Type_commit(&elementType);
201
202 /*commit other type */
203 offsets[0] = 0;
204 oldtypes[0] = MPI_INT;
205 blockcounts[0] = 15;
206
207 offsets[1] = 15 * int_ext;
208 oldtypes[1] = MPI_SHORT;
209 blockcounts[1] = 4;
210
211 blockcounts[2] = 1;
212 oldtypes[2] = MPI_UNSIGNED;
213 offsets[2] = offsets[1] + 4 * short_ext;
214
215 blockcounts[3] = 1024 + 1024 + 1024;
216 oldtypes[3] = MPI_CHAR;
217 offsets[3] = offsets[2] + 1 * uint_ext;
218
219 /* Now define structured type and commit it */
220 MPI_Type_create_struct(4, blockcounts, offsets, oldtypes, &otherType);
221 MPI_Type_commit(&otherType);
222
223 if (MPI_dataset->myid == 0) {
224 /*fill other layout structure */
225 other.n_columns = layout->n_columns;
226 other.n_parameters = layout->n_parameters;
227 other.n_associates = layout->n_associates;
228 other.n_arrays = layout->n_arrays;
229 other.version = layout->version;
230 other.layout_offset = SDDS_dataset->pagecount_offset[0];
231 other.layout_written = layout->layout_written;
232 other.disconnected = layout->disconnected;
233 other.gzipFile = layout->gzipFile;
234 other.lzmaFile = layout->lzmaFile;
235 other.popenUsed = layout->popenUsed;
236 other.depth = layout->depth;
237 other.data_command_seen = layout->data_command_seen;
238 other.commentFlags = layout->commentFlags;
239 other.byteOrderDeclared = layout->byteOrderDeclared;
240 other.mode = layout->data_mode.mode;
241 other.lines_per_row = layout->data_mode.lines_per_row;
242 other.no_row_counts = layout->data_mode.no_row_counts;
243 other.fixed_row_count = layout->data_mode.fixed_row_count;
244 other.column_memory_mode = layout->data_mode.column_memory_mode;
245 other.fsync_data = layout->data_mode.fsync_data;
246 other.additional_header_lines = layout->data_mode.additional_header_lines;
247 other.description_len = other.contents_len = other.filename_len = 0;
248 other.description[0] = other.contents[0] = other.filename[0] = '\0';
249 other.swapByteOrder = SDDS_dataset->swapByteOrder;
250
251 if (layout->description) {
252 other.description_len = strlen(layout->description);
253 sprintf(other.description, "%s", layout->description);
254 }
255 if (layout->contents) {
256 other.contents_len = strlen(layout->contents);
257 sprintf(other.contents, "%s", layout->contents);
258 }
259 if (layout->filename) {
260 other.filename_len = strlen(layout->filename);
261 sprintf(other.filename, "%s", layout->filename);
262 }
263 } else {
264 /* if (!SDDS_ZeroMemory((void *)SDDS_dataset, sizeof(SDDS_DATASET))) {
265 SDDS_SetError("Unable to zero memory for sdds dataset.");
266 return 0;
267 }
268 SDDS_dataset->MPI_dataset = MPI_dataset; */
269 }
270 /* broadcaset the layout other to other processors */
271 MPI_Bcast(&other, 1, otherType, 0, MPI_dataset->comm);
272 MPI_Type_free(&otherType);
273 if (other.n_columns)
274 column = malloc(sizeof(*column) * other.n_columns);
275 if (other.n_parameters)
276 parameter = malloc(sizeof(*parameter) * other.n_parameters);
277 if (other.n_arrays)
278 array = malloc(sizeof(*array) * other.n_arrays);
279 if (other.n_associates)
280 associate = malloc(sizeof(*associate) * other.n_associates);
281 if (MPI_dataset->myid == 0) {
282 /*fill elements */
283 for (i = 0; i < other.n_columns; i++) {
284 if (!SDDS_ZeroMemory(&column[i], sizeof(ELEMENT_DEF))) {
285 SDDS_SetError("Unable to zero memory for columns(SDDS_MPI_InitializeInput)");
286 return 0;
287 }
288 column[i].type = layout->column_definition[i].type;
289 column[i].field_length = layout->column_definition[i].field_length;
290 column[i].definition_mode = layout->column_definition[i].definition_mode;
291 column[i].memory_number = layout->column_definition[i].memory_number;
292 column[i].pointer_number = layout->column_definition[i].pointer_number;
293 if (layout->column_definition[i].name) {
294 column[i].name_len = strlen(layout->column_definition[i].name);
295 sprintf(column[i].name, "%s", layout->column_definition[i].name);
296 }
297 if (layout->column_definition[i].symbol) {
298 column[i].symbol_len = strlen(layout->column_definition[i].symbol);
299 sprintf(column[i].symbol, "%s", layout->column_definition[i].symbol);
300 }
301 if (layout->column_definition[i].units) {
302 column[i].units_len = strlen(layout->column_definition[i].units);
303 sprintf(column[i].units, "%s", layout->column_definition[i].units);
304 }
305 if (layout->column_definition[i].description) {
306 column[i].description_len = strlen(layout->column_definition[i].description);
307 sprintf(column[i].description, "%s", layout->column_definition[i].description);
308 }
309 if (layout->column_definition[i].format_string) {
310 column[i].format_string_len = strlen(layout->column_definition[i].format_string);
311 sprintf(column[i].format_string, "%s", layout->column_definition[i].format_string);
312 }
313 }
314 for (i = 0; i < other.n_parameters; i++) {
315 if (!SDDS_ZeroMemory(&parameter[i], sizeof(ELEMENT_DEF))) {
316 SDDS_SetError("Unable to zero memory for parameters(SDDS_MPI_InitializeInput)");
317 return 0;
318 }
319 parameter[i].type = layout->parameter_definition[i].type;
320 parameter[i].definition_mode = layout->parameter_definition[i].definition_mode;
321 parameter[i].memory_number = layout->parameter_definition[i].memory_number;
322 if (layout->parameter_definition[i].name) {
323 parameter[i].name_len = strlen(layout->parameter_definition[i].name);
324 sprintf(parameter[i].name, "%s", layout->parameter_definition[i].name);
325 }
326 if (layout->parameter_definition[i].symbol) {
327 parameter[i].symbol_len = strlen(layout->parameter_definition[i].symbol);
328 sprintf(parameter[i].symbol, "%s", layout->parameter_definition[i].symbol);
329 }
330 if (layout->parameter_definition[i].units) {
331 parameter[i].units_len = strlen(layout->parameter_definition[i].units);
332 sprintf(parameter[i].units, "%s", layout->parameter_definition[i].units);
333 }
334 if (layout->parameter_definition[i].description) {
335 parameter[i].description_len = strlen(layout->parameter_definition[i].description);
336 sprintf(parameter[i].description, "%s", layout->parameter_definition[i].description);
337 }
338 if (layout->parameter_definition[i].format_string) {
339 parameter[i].format_string_len = strlen(layout->parameter_definition[i].format_string);
340 sprintf(parameter[i].format_string, "%s", layout->parameter_definition[i].format_string);
341 }
342 if (layout->parameter_definition[i].fixed_value) {
343 parameter[i].fixed_value_len = strlen(layout->parameter_definition[i].fixed_value);
344 sprintf(parameter[i].fixed_value, "%s", layout->parameter_definition[i].fixed_value);
345 } else
346 parameter[i].fixed_value_len = -1;
347 }
348 for (i = 0; i < other.n_arrays; i++) {
349 if (!SDDS_ZeroMemory(&array[i], sizeof(ELEMENT_DEF))) {
350 SDDS_SetError("Unable to zero memory for arrays(SDDS_MPI_InitializeInput)");
351 return 0;
352 }
353 array[i].type = layout->array_definition[i].type;
354 array[i].field_length = layout->array_definition[i].field_length;
355 array[i].dimensions = layout->array_definition[i].dimensions;
356 if (layout->array_definition[i].name) {
357 array[i].name_len = strlen(layout->array_definition[i].name);
358 sprintf(array[i].name, "%s", layout->array_definition[i].name);
359 }
360 if (layout->array_definition[i].symbol) {
361 array[i].symbol_len = strlen(layout->array_definition[i].symbol);
362 sprintf(array[i].symbol, "%s", layout->array_definition[i].symbol);
363 }
364 if (layout->array_definition[i].units) {
365 array[i].units_len = strlen(layout->array_definition[i].units);
366 sprintf(array[i].units, "%s", layout->array_definition[i].units);
367 }
368 if (layout->array_definition[i].description) {
369 array[i].description_len = strlen(layout->array_definition[i].description);
370 sprintf(array[i].description, "%s", layout->array_definition[i].description);
371 }
372 if (layout->array_definition[i].format_string) {
373 array[i].format_string_len = strlen(layout->array_definition[i].format_string);
374 sprintf(array[i].format_string, "%s", layout->array_definition[i].format_string);
375 }
376 if (layout->array_definition[i].group_name) {
377 array[i].group_name_len = strlen(layout->array_definition[i].group_name);
378 sprintf(array[i].group_name, "%s", layout->array_definition[i].group_name);
379 }
380 }
381 for (i = 0; i < other.n_associates; i++) {
382 if (!SDDS_ZeroMemory(&associate[i], sizeof(ASSOCIATE_DEF))) {
383 SDDS_SetError("Unable to zero memory for associates(SDDS_MPI_InitializeInput)");
384 return 0;
385 }
386 associate[i].sdds = layout->associate_definition[i].sdds;
387
388 if (layout->associate_definition[i].name) {
389 associate[i].name_len = strlen(layout->associate_definition[i].name);
390 sprintf(associate[i].name, "%s", layout->associate_definition[i].name);
391 }
392 if (layout->associate_definition[i].filename) {
393 associate[i].filename_len = strlen(layout->associate_definition[i].filename);
394 sprintf(associate[i].filename, "%s", layout->associate_definition[i].filename);
395 }
396 if (layout->associate_definition[i].path) {
397 associate[i].path_len = strlen(layout->associate_definition[i].path);
398 sprintf(associate[i].path, "%s", layout->associate_definition[i].path);
399 }
400 if (layout->associate_definition[i].description) {
401 associate[i].description_len = strlen(layout->associate_definition[i].description);
402 sprintf(associate[i].description, "%s", layout->associate_definition[i].description);
403 }
404 if (layout->associate_definition[i].contents) {
405 associate[i].contents_len = strlen(layout->associate_definition[i].contents);
406 sprintf(associate[i].contents, "%s", layout->associate_definition[i].contents);
407 }
408 }
409 } else {
410 SDDS_dataset->page_number = SDDS_dataset->page_started = 0;
411 layout->popenUsed = other.popenUsed;
412 layout->gzipFile = other.gzipFile;
413 layout->lzmaFile = other.lzmaFile;
414 layout->depth = other.depth;
415 layout->data_command_seen = other.data_command_seen;
416 layout->commentFlags = other.commentFlags;
417 layout->disconnected = other.disconnected;
418 layout->layout_written = other.layout_written;
419 if (other.filename_len)
420 SDDS_CopyString(&layout->filename, other.filename);
421 layout->version = other.version;
422 layout->byteOrderDeclared = other.byteOrderDeclared;
423 layout->data_mode.mode = other.mode;
424 layout->data_mode.lines_per_row = other.lines_per_row;
425 layout->data_mode.no_row_counts = other.no_row_counts;
426 layout->data_mode.fixed_row_count = other.fixed_row_count;
427 layout->data_mode.fsync_data = other.fsync_data;
428 layout->data_mode.column_memory_mode = other.column_memory_mode;
429 layout->data_mode.additional_header_lines = other.additional_header_lines;
430 if (other.description_len)
431 SDDS_CopyString(&layout->description, other.description);
432 if (other.contents_len)
433 SDDS_CopyString(&layout->contents, other.contents);
434 SDDS_dataset->swapByteOrder = other.swapByteOrder;
435 }
436 if (other.n_columns)
437 MPI_Bcast(column, other.n_columns, elementType, 0, MPI_dataset->comm);
438 if (other.n_parameters)
439 MPI_Bcast(parameter, other.n_parameters, elementType, 0, MPI_dataset->comm);
440 if (other.n_arrays)
441 MPI_Bcast(array, other.n_arrays, elementType, 0, MPI_dataset->comm);
442 MPI_Type_free(&elementType);
443 if (other.n_associates) {
444 /* create and commit associate type */
445 offsets[0] = 0;
446 oldtypes[0] = MPI_INT;
447 blockcounts[0] = 6;
448
449 offsets[1] = 6 * int_ext;
450 oldtypes[1] = MPI_CHAR;
451 blockcounts[1] = 256 + 256 + 1024 + 1024 + 1024;
452
453 MPI_Type_create_struct(2, blockcounts, offsets, oldtypes, &associateType);
454 MPI_Type_commit(&associateType);
455 MPI_Bcast(associate, other.n_associates, associateType, 0, MPI_dataset->comm);
456 MPI_Type_free(&associateType);
457 }
458 if (MPI_dataset->myid) {
459 for (i = 0; i < other.n_columns; i++) {
460 symbol = units = description = format_string = NULL;
461 if (column[i].units_len)
462 SDDS_CopyString(&units, column[i].units);
463 if (column[i].description_len)
464 SDDS_CopyString(&description, column[i].description);
465 if (column[i].format_string_len)
466 SDDS_CopyString(&format_string, column[i].format_string);
467 if (column[i].symbol_len)
468 SDDS_CopyString(&symbol, column[i].symbol);
469 if (SDDS_DefineColumn(SDDS_dataset, column[i].name, symbol, units, description, format_string, column[i].type, column[i].field_length) < 0) {
470 SDDS_SetError("Unable to define column (SDDS_MPI_BroadcastLayout)");
471 return (0);
472 }
473 if (units)
474 free(units);
475 if (description)
476 free(description);
477 if (symbol)
478 free(symbol);
479 if (format_string)
480 free(format_string);
481 }
482 for (i = 0; i < other.n_parameters; i++) {
483 symbol = units = description = format_string = fixed_value = NULL;
484 if (parameter[i].units_len)
485 SDDS_CopyString(&units, parameter[i].units);
486 if (parameter[i].description_len)
487 SDDS_CopyString(&description, parameter[i].description);
488 if (parameter[i].format_string_len)
489 SDDS_CopyString(&format_string, parameter[i].format_string);
490 if (parameter[i].symbol_len)
491 SDDS_CopyString(&symbol, parameter[i].symbol);
492 if (parameter[i].fixed_value_len >= 0)
493 SDDS_CopyString(&fixed_value, parameter[i].fixed_value);
494 if (SDDS_DefineParameter(SDDS_dataset, parameter[i].name, symbol, units, description, format_string, parameter[i].type, fixed_value) < 0) {
495 SDDS_SetError("Unable to define parameter (SDDS_MPI_BroadcastLayout)");
496 return (0);
497 }
498 if (units)
499 free(units);
500 if (description)
501 free(description);
502 if (symbol)
503 free(symbol);
504 if (format_string)
505 free(format_string);
506 if (fixed_value)
507 free(fixed_value);
508 }
509 for (i = 0; i < other.n_arrays; i++) {
510 if (SDDS_DefineArray(SDDS_dataset, array[i].name, array[i].symbol, array[i].units, array[i].description, array[i].format_string, array[i].type, array[i].field_length, array[i].dimensions, array[i].group_name) < 0) {
511 SDDS_SetError("Unable to define array (SDDS_BroadcastLayout)");
512 return (0);
513 }
514 }
515 for (i = 0; i < other.n_associates; i++) {
516 filename = path = description = contents = NULL;
517 if (associate[i].filename_len)
518 SDDS_CopyString(&filename, associate[i].filename);
519 if (associate[i].path_len)
520 SDDS_CopyString(&path, associate[i].path);
521 if (associate[i].description_len)
522 SDDS_CopyString(&description, associate[i].description);
523 if (associate[i].contents_len)
524 SDDS_CopyString(&contents, associate[i].contents);
525 if (SDDS_DefineAssociate(SDDS_dataset, associate[i].name, filename, path, description, contents, associate[i].sdds) < 0) {
526 SDDS_SetError("Unable to define associate (SDDS_MPI_BroadcastLayout)");
527 return (0);
528 }
529 if (description)
530 free(description);
531 if (contents)
532 free(contents);
533 if (path)
534 free(path);
535 if (filename)
536 free(filename);
537 }
538 if (!SDDS_SaveLayout(SDDS_dataset)) {
539 SDDS_SetError("Unable to save layout (SDDS_BroadcastLayout)");
540 return (0);
541 }
542 }
543 if (column)
544 free(column);
545 if (array)
546 free(array);
547 if (parameter)
548 free(parameter);
549 if (associate)
550 free(associate);
551 column = array = parameter = NULL;
552 associate = NULL;
553 MPI_dataset->file_offset = other.layout_offset;
554 return 1;
555}
int32_t SDDS_SaveLayout(SDDS_DATASET *SDDS_dataset)
Definition SDDS_copy.c:615
int32_t SDDS_DefineAssociate(SDDS_DATASET *SDDS_dataset, const char *name, const char *filename, const char *path, const char *description, const char *contents, int32_t sdds)
Defines an associate for the SDDS dataset.
int32_t SDDS_DefineArray(SDDS_DATASET *SDDS_dataset, const char *name, const char *symbol, const char *units, const char *description, const char *format_string, int32_t type, int32_t field_length, int32_t dimensions, const char *group_name)
Defines a data array within the SDDS dataset.
int32_t SDDS_DefineColumn(SDDS_DATASET *SDDS_dataset, const char *name, const char *symbol, const char *units, const char *description, const char *format_string, int32_t type, int32_t field_length)
Defines a data column within the SDDS dataset.
int32_t SDDS_DefineParameter(SDDS_DATASET *SDDS_dataset, const char *name, const char *symbol, const char *units, const char *description, const char *format_string, int32_t type, char *fixed_value)
Defines a data parameter with a fixed string value.
Structure defining an associate with various attributes.
Structure defining an element with various attributes.
Structure defining additional layout information.

◆ SDDS_MPI_InitializeInput()

int32_t SDDS_MPI_InitializeInput ( SDDS_DATASET * SDDS_dataset,
char * filename )

Initializes an SDDS dataset for input using MPI.

This function initializes the provided SDDS_DATASET structure for reading data from a specified file. It handles various file formats, including plain, gzip, and LZMA-compressed files. The function sets up necessary layout information, allocates memory for columns and related structures, and prepares the dataset for parallel I/O operations using MPI.

Parameters
[in,out]SDDS_datasetPointer to the SDDS_DATASET structure to be initialized.
[in]filenameThe name of the file to be opened for reading. If filename is NULL, the function will attempt to read from the standard input (stdin).
Returns
  • 1 on successful initialization.
  • 0 if an error occurs during initialization, such as memory allocation failures, file opening issues, or layout reading errors.

The function performs the following steps:

  • Validates the SDDS_DATASET structure.
  • Initializes layout-related flags and variables.
  • Copies the filename into the dataset's layout if provided.
  • Opens the specified file, handling different compression formats based on the file extension (.gz, .lzma, .xz).
  • Reads the layout information from the file, supporting gzip and LZMA compression.
  • Allocates memory for column flags and orders if columns are defined.
  • Sets the dataset mode to read-only and initializes file offsets for sequential reading.
  • Closes the file pointer after reading the layout.
  • Broadcasts the layout information to all MPI processes if MASTER_READTITLE_ONLY is defined.
  • Opens the MPI file for parallel reading and retrieves the file size and column offsets.
Note
  • The function assumes that only one page of data is present in the input file.
  • If MASTER_READTITLE_ONLY is defined, only the root process reads the title.
See also
SDDS_CheckDataset, SDDS_SetError, SDDS_CopyString, SDDS_ReadLayout, SDDS_GZipReadLayout, SDDS_LZMAReadLayout, SDDS_SaveLayout, SDDS_MPI_BroadcastLayout, SDDS_MPI_File_Open, SDDS_MPI_Get_Column_Size

Definition at line 601 of file SDDSmpi_input.c.

601 {
602 /* startTime = MPI_Wtime(); */
603 MPI_DATASET *MPI_dataset = SDDS_dataset->MPI_dataset;
604
605#if defined(MASTER_READTITLE_ONLY)
606 if (MPI_dataset->myid == 0)
607#endif
608 {
609 /* char *ptr, *datafile, *headerfile; */
610 static char s[SDDS_MAXLINE];
611#if defined(zLib)
612 char *extension;
613#endif
614 if (sizeof(gzFile) != sizeof(void *)) {
615 SDDS_SetError("gzFile is not the same size as void *, possible corruption of the SDDS_LAYOUT structure");
616 return (0);
617 }
618 if (!SDDS_CheckDataset(SDDS_dataset, "SDDS_InitializeInput"))
619 return (0);
620 SDDS_dataset->layout.gzipFile = SDDS_dataset->layout.lzmaFile = SDDS_dataset->layout.disconnected = SDDS_dataset->layout.popenUsed = 0;
621 SDDS_dataset->layout.depth = SDDS_dataset->layout.data_command_seen = SDDS_dataset->layout.commentFlags = 0;
622 if (!filename)
623 SDDS_dataset->layout.filename = NULL;
624 else if (!SDDS_CopyString(&SDDS_dataset->layout.filename, filename)) {
625 sprintf(s, "Memory allocation failure initializing file \"%s\" (SDDS_InitializeInput)", filename);
626 SDDS_SetError(s);
627 return (0);
628 }
629 if (!filename) {
630#if defined(_WIN32)
631 if (_setmode(_fileno(stdin), _O_BINARY) == -1) {
632 sprintf(s, "unable to set stdin to binary mode");
633 SDDS_SetError(s);
634 return 0;
635 }
636#endif
637 SDDS_dataset->layout.fp = stdin;
638 } else {
639#if defined(zLib)
640 if (!(extension = strrchr(filename, '.')) || strcmp(extension, ".gz") != 0) {
641#endif
642 if ((extension = strrchr(filename, '.')) && ((strcmp(extension, ".lzma") == 0) || (strcmp(extension, ".xz") == 0))) {
643 SDDS_dataset->layout.lzmaFile = 1;
644 if (!(SDDS_dataset->layout.lzmafp = UnpackLZMAOpen(filename))) {
645 sprintf(s, "Unable to open file \"%s\" for reading (SDDS_InitializeInput)", filename);
646 SDDS_SetError(s);
647 return (0);
648 }
649 SDDS_dataset->layout.fp = SDDS_dataset->layout.lzmafp->fp;
650 } else {
651 if (!(SDDS_dataset->layout.fp = UnpackFopen(filename, UNPACK_REQUIRE_SDDS | UNPACK_USE_PIPE, &SDDS_dataset->layout.popenUsed, NULL))) {
652 sprintf(s, "Unable to open file \"%s\" for reading (SDDS_InitializeInput)", filename);
653 SDDS_SetError(s);
654 return (0);
655 }
656 }
657#if defined(zLib)
658 } else {
659 SDDS_dataset->layout.gzipFile = 1;
660 if (!(SDDS_dataset->layout.gzfp = UnpackGZipOpen(filename))) {
661 sprintf(s, "Unable to open file \"%s\" for reading (SDDS_InitializeInput)", filename);
662 SDDS_SetError(s);
663 return (0);
664 }
665 }
666#endif
667 }
668 SDDS_dataset->page_number = SDDS_dataset->page_started = 0;
669 SDDS_dataset->file_had_data = 0;
670 SDDS_DeferSavingLayout(SDDS_dataset, 1);
671#if defined(zLib)
672 if (SDDS_dataset->layout.gzipFile) {
673 if (!SDDS_GZipReadLayout(SDDS_dataset, SDDS_dataset->layout.gzfp))
674 return (0);
675 } else {
676#endif
677 if (SDDS_dataset->layout.lzmaFile) {
678 if (!SDDS_LZMAReadLayout(SDDS_dataset, SDDS_dataset->layout.lzmafp))
679 return (0);
680 } else {
681 if (!SDDS_ReadLayout(SDDS_dataset, SDDS_dataset->layout.fp))
682 return (0);
683 }
684#if defined(zLib)
685 }
686#endif
687 SDDS_dataset->layout.layout_written = 0;
688 SDDS_DeferSavingLayout(SDDS_dataset, 0);
689 if (!SDDS_SaveLayout(SDDS_dataset))
690 return 0;
691 if (SDDS_dataset->layout.n_columns &&
692 ((!(SDDS_dataset->column_flag = (int32_t *)SDDS_Malloc(sizeof(int32_t) * SDDS_dataset->layout.n_columns)) ||
693 !(SDDS_dataset->column_order = (int32_t *)SDDS_Malloc(sizeof(int32_t) * SDDS_dataset->layout.n_columns))) ||
694 (!SDDS_SetMemory(SDDS_dataset->column_flag, SDDS_dataset->layout.n_columns, SDDS_LONG, (int32_t)1, (int32_t)0) || !SDDS_SetMemory(SDDS_dataset->column_order, SDDS_dataset->layout.n_columns, SDDS_LONG, (int32_t)0, (int32_t)1)))) {
695 SDDS_SetError("Unable to initialize input--memory allocation failure (SDDS_InitializeInput)");
696 return (0);
697 }
698 SDDS_dataset->mode = SDDS_READMODE; /*reading */
699 SDDS_dataset->pagecount_offset = NULL;
700 if (!SDDS_dataset->layout.gzipFile && !SDDS_dataset->layout.lzmaFile && !SDDS_dataset->layout.popenUsed && SDDS_dataset->layout.filename) {
701 /* Data is not:
702 1. from a gzip file
703 2. from a file that is being internally decompressed by a command executed with popen()
704 3. from a pipe set up externally (e.g., -pipe=in on commandline)
705 */
706 SDDS_dataset->pages_read = 0;
707 SDDS_dataset->pagecount_offset = malloc(sizeof(*SDDS_dataset->pagecount_offset));
708 SDDS_dataset->pagecount_offset[0] = ftell(SDDS_dataset->layout.fp);
709 fseek(SDDS_dataset->layout.fp, 0, 2); /*point to the end of the file */
710 SDDS_dataset->endOfFile_offset = ftell(SDDS_dataset->layout.fp);
711 fseek(SDDS_dataset->layout.fp, SDDS_dataset->pagecount_offset[0], 0);
712 /*point to the beginning of the first page */
713 }
714 if (SDDS_dataset->layout.fp)
715 fclose(SDDS_dataset->layout.fp);
716 /* SDDS_dataset->MPI_dataset = MPI_dataset; */
717 }
718
719#if defined(MASTER_READTITLE_ONLY)
720 if (!SDDS_MPI_BroadcastLayout(SDDS_dataset))
721 return 0;
722#else
723 MPI_dataset->file_offset = SDDS_dataset->pagecount_offset[0];
724#endif
725
726 if (!SDDS_MPI_File_Open(MPI_dataset, filename, SDDS_MPI_READ_ONLY)) {
727 SDDS_SetError("(SDDS_MPI_File_Open)Unablle to open file for reading.");
728 return 0;
729 }
730
731 MPI_File_get_size(MPI_dataset->MPI_file, &(MPI_dataset->file_size));
732 if (MPI_dataset->file_offset >= MPI_dataset->file_size) {
733 SDDS_SetError("No data contained in the input file.");
734 return 0;
735 }
736 MPI_dataset->column_offset = SDDS_MPI_Get_Column_Size(SDDS_dataset);
737 SDDS_dataset->parallel_io = 1;
738 return 1;
739}
MPI_Offset SDDS_MPI_Get_Column_Size(SDDS_DATASET *SDDS_dataset)
Get the total size of all columns in an SDDS dataset.
void SDDS_DeferSavingLayout(SDDS_DATASET *SDDS_dataset, int32_t mode)
Definition SDDS_copy.c:603
int32_t SDDS_ReadLayout(SDDS_DATASET *SDDS_dataset, FILE *fp)
Definition SDDS_input.c:517
int32_t SDDS_LZMAReadLayout(SDDS_DATASET *SDDS_dataset, struct lzmafile *lzmafp)
Definition SDDS_input.c:680
int32_t SDDS_SetMemory(void *mem, int64_t n_elements, int32_t data_type,...)
Initializes a memory block with a sequence of values based on a specified data type.
int32_t SDDS_CheckDataset(SDDS_DATASET *SDDS_dataset, const char *caller)
Validates the SDDS dataset pointer.
Definition SDDS_utils.c:552
void * SDDS_Malloc(size_t size)
Allocates memory of a specified size.
Definition SDDS_utils.c:639
int32_t SDDS_MPI_File_Open(MPI_DATASET *MPI_dataset, char *filename, unsigned long flags)
Opens an MPI file with the specified flags.

◆ SDDS_MPI_InitializeInputFromSearchPath()

int32_t SDDS_MPI_InitializeInputFromSearchPath ( SDDS_DATASET * SDDS_dataset,
char * file )

Initializes an SDDS dataset for input by searching the provided search path using MPI.

This function searches for the specified file within predefined search paths and initializes the SDDS_DATASET structure for reading using MPI. It leverages SDDS_MPI_InitializeInput to perform the actual initialization once the file is found.

Parameters
[in,out]SDDS_datasetPointer to the SDDS_DATASET structure to be initialized.
[in]fileThe name of the file to search for and initialize. The function searches for this file in the configured search paths.
Returns
  • 1 on successful initialization.
  • 0 if the file is not found in the search paths or if initialization fails.

The function performs the following steps:

  • Searches for the specified file in the predefined search paths using findFileInSearchPath.
  • If the file is found, it calls SDDS_MPI_InitializeInput to initialize the dataset.
  • Frees the allocated memory for the filename after initialization.
See also
SDDS_MPI_InitializeInput, findFileInSearchPath, SDDS_SetError

Definition at line 767 of file SDDSmpi_input.c.

767 {
768 char *filename;
769 int32_t value;
770 if (!(filename = findFileInSearchPath(file))) {
771 SDDS_SetError("file does not exist in searchpath (InitializeInputFromSearchPath)");
772 return 0;
773 }
774 value = SDDS_MPI_InitializeInput(SDDS_dataset, filename);
775 free(filename);
776 return value;
777}
int32_t SDDS_MPI_InitializeInput(SDDS_DATASET *SDDS_dataset, char *filename)
Initializes an SDDS dataset for input using MPI.

◆ SDDS_MPI_ReadPage()

int32_t SDDS_MPI_ReadPage ( SDDS_DATASET * SDDS_dataset)

Reads a page from an SDDS dataset using MPI.

This function checks the validity of the provided SDDS dataset and ensures that the dataset is connected and in binary mode before attempting to read a page. If the dataset is in ASCII mode or disconnected, appropriate errors are set and the function returns 0. On successful reading of a binary page, the function returns 1.

Parameters
SDDS_datasetPointer to the SDDS_DATASET structure representing the dataset.
Returns
  • 1 on successful page read.
  • 0 if the dataset is invalid, disconnected, in ASCII mode, or an unrecognized data mode.
See also
SDDS_CheckDataset, SDDS_SetError, SDDS_MPI_ReadBinaryPage

Definition at line 125 of file SDDSmpi_input.c.

125 {
126 if (!SDDS_CheckDataset(SDDS_dataset, "SDDS_ReadPageSparse"))
127 return (0);
128 if (SDDS_dataset->layout.disconnected) {
129 SDDS_SetError("Can't read page--file is disconnected (SDDS_ReadPageSparse)");
130 return 0;
131 }
132 if (SDDS_dataset->original_layout.data_mode.mode == SDDS_ASCII) {
133 SDDS_SetError("Unable to read ascii file with SDDS_MPI.");
134 return 0;
135 } else if (SDDS_dataset->original_layout.data_mode.mode == SDDS_BINARY) {
136 return SDDS_MPI_ReadBinaryPage(SDDS_dataset);
137 } else {
138 SDDS_SetError("Unable to read page--unrecognized data mode (SDDS_ReadPageSparse)");
139 return (0);
140 }
141}
int32_t SDDS_MPI_ReadBinaryPage(SDDS_DATASET *SDDS_dataset)
Reads a binary page from an SDDS dataset using MPI parallel I/O.