File: bpWriter.c

package info (click to toggle)
adios2 2.10.2%2Bdfsg1-2
  • links: PTS, VCS
  • area: main
  • in suites: trixie
  • size: 33,764 kB
  • sloc: cpp: 175,964; ansic: 160,510; f90: 14,630; yacc: 12,668; python: 7,275; perl: 7,126; sh: 2,825; lisp: 1,106; xml: 1,049; makefile: 579; lex: 557
file content (115 lines) | stat: -rw-r--r-- 2,594 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
/*
 * Distributed under the OSI-approved Apache License, Version 2.0.  See
 * accompanying file Copyright.txt for details.
 *
 * bpWriter.c : C bindings version of bpWriter.cpp
 *
 *  Created on: Aug 8, 2017
 *      Author: William F Godoy godoywf@ornl.gov
 */
#include <stdio.h>  // printf
#include <stdlib.h> // malloc, free, exit

#include <adios2_c.h>

#if ADIOS2_USE_MPI
#include <mpi.h>
#endif

void check_error(const int error)
{
    if (error)
    {
        printf("ERROR: %d\n", error);
        exit(error);
    }
}

void check_handler(const void *handler, const char *message)
{
    if (handler == NULL)
    {
        printf("ERROR: invalid %s handler \n", message);
        exit(EXIT_FAILURE);
    }
}

int main(int argc, char *argv[])
{
    int rank, size;

#if ADIOS2_USE_MPI
    int provided;

    // MPI_THREAD_MULTIPLE is only required if you enable the SST MPI_DP
    MPI_Init_thread(&argc, &argv, MPI_THREAD_MULTIPLE, &provided);
    MPI_Comm_rank(MPI_COMM_WORLD, &rank);
    MPI_Comm_size(MPI_COMM_WORLD, &size);
#else
    rank = 0;
    size = 1;
#endif

    adios2_error errio;
    // application input, data in heap
    const size_t Nx = 10;
    float *myFloats;
    myFloats = malloc(sizeof(float) * Nx);

    unsigned int i;
    for (i = 0; i < Nx; ++i)
    {
        myFloats[i] = (float)i;
    }

#if ADIOS2_USE_MPI
    adios2_adios *adios = adios2_init_mpi(MPI_COMM_WORLD);
#else
    adios2_adios *adios = adios2_init_serial();
#endif

    adios2_step_status err;
    check_handler(adios, "adios");

    adios2_io *io = adios2_declare_io(adios, "BPFile_Write");
    check_handler(io, "io");

    // dims are allocated in the stack
    size_t shape[1];
    shape[0] = (size_t)size * Nx;

    size_t start[1];
    start[0] = (size_t)rank * Nx;

    size_t count[1];
    count[0] = Nx;

    adios2_variable *variable = adios2_define_variable(io, "bpFloats", adios2_type_float, 1, shape,
                                                       start, count, adios2_constant_dims_true);
    check_handler(variable, "variable");

    adios2_engine *engine = adios2_open(io, "myVector_c.bp", adios2_mode_write);
    check_handler(engine, "engine");

    adios2_begin_step(engine, adios2_step_mode_append, 0.0f, &err);

    errio = adios2_put(engine, variable, myFloats, adios2_mode_deferred);
    check_error(errio);

    adios2_end_step(engine);

    errio = adios2_close(engine);
    check_error(errio);

    // deallocate adios
    errio = adios2_finalize(adios);
    check_error(errio);

    free(myFloats);

#if ADIOS2_USE_MPI
    MPI_Finalize();
#endif

    return 0;
}