File: win_shared_zerobyte.c

package info (click to toggle)
mpich 4.0.2-3
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 423,384 kB
  • sloc: ansic: 1,088,434; cpp: 71,364; javascript: 40,763; f90: 22,829; sh: 17,463; perl: 14,773; xml: 14,418; python: 10,265; makefile: 9,246; fortran: 8,008; java: 4,355; asm: 324; ruby: 176; lisp: 19; php: 8; sed: 4
file content (117 lines) | stat: -rw-r--r-- 3,185 bytes parent folder | download | duplicates (4)
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
116
117
/*
 * Copyright (C) by Argonne National Laboratory
 *     See COPYRIGHT in top-level directory
 */

#include <stdio.h>
#include <stdlib.h>
#include <assert.h>
#include <mpi.h>
#include "mpitest.h"

const int verbose = 0;

int main(int argc, char **argv)
{
    int i, rank, nproc;
    int shm_rank, shm_nproc, last_shm_rank_w;
    MPI_Aint size;
    int errors = 0;
    int **bases = NULL, *abs_base, *my_base;
    int disp_unit;
    MPI_Win shm_win;
    MPI_Comm shm_comm;
    int ELEM_PER_PROC = 0;

    MTest_Init(&argc, &argv);

    MPI_Comm_rank(MPI_COMM_WORLD, &rank);
    MPI_Comm_size(MPI_COMM_WORLD, &nproc);

    MPI_Comm_split_type(MPI_COMM_WORLD, MPI_COMM_TYPE_SHARED, rank, MPI_INFO_NULL, &shm_comm);

    MPI_Comm_rank(shm_comm, &shm_rank);
    MPI_Comm_size(shm_comm, &shm_nproc);

    /* Platform does not support shared memory, just return. */
    if (shm_nproc < 2) {
        goto exit;
    }

    bases = calloc(shm_nproc, sizeof(int *));

    if (shm_rank == 0 || shm_rank == shm_nproc - 1) {
        ELEM_PER_PROC = 16;
    }

    /* Obtain the world rank of last process in my shm_comm. */
    if (shm_rank == shm_nproc - 1)
        last_shm_rank_w = rank;
    MPI_Bcast(&last_shm_rank_w, 1, MPI_INT, shm_nproc - 1, shm_comm);
    if (verbose)
        printf("%d -- my rank %d, last shm rank %d in world \n", shm_rank, rank, last_shm_rank_w);

    /* Allocate ELEM_PER_PROC integers for each process */
    MPI_Win_allocate_shared(sizeof(int) * ELEM_PER_PROC, sizeof(int), MPI_INFO_NULL,
                            shm_comm, &my_base, &shm_win);

    /* Locate absolute base */
    MPI_Win_shared_query(shm_win, MPI_PROC_NULL, &size, &disp_unit, &abs_base);

    if (verbose)
        printf("%d -- allocate shared: my_base = %p, absolute base = %p\n", shm_rank, my_base,
               abs_base);

    for (i = 0; i < shm_nproc; i++) {
        MPI_Win_shared_query(shm_win, i, &size, &disp_unit, &bases[i]);
        if (verbose)
            printf("%d --    shared query: base[%d]=%p, size %ld, unit %d\n",
                   shm_rank, i, bases[i], size, disp_unit);
    }

    MPI_Win_lock_all(MPI_MODE_NOCHECK, shm_win);

    /* Reset data by using my world rank */
    for (i = 0; i < ELEM_PER_PROC; i++) {
        my_base[i] = rank * -1;
    }

    MPI_Win_sync(shm_win);
    MPI_Barrier(shm_comm);
    MPI_Win_sync(shm_win);

    /* Read and verify everyone's data */
    if (shm_rank == 0) {
        bases[0][0] = 1;
    }

    MPI_Win_sync(shm_win);
    MPI_Barrier(shm_comm);
    MPI_Win_sync(shm_win);

    if (bases[0][0] != 1) {
        errors++;
        printf("%d -- my rank %d, got %d at rank %d index %d, expected %d\n", shm_rank, rank,
               bases[0][0], 0, 0, 1);
    }

    if (bases[shm_nproc - 1][0] != last_shm_rank_w * -1) {
        errors++;
        printf("%d -- my rank %d, got %d at rank %d index %d, expected %d\n", shm_rank, rank,
               bases[shm_nproc - 1][0], shm_nproc - 1, 0, last_shm_rank_w * -1);
    }

    MPI_Win_unlock_all(shm_win);
    MPI_Win_free(&shm_win);

  exit:

    MPI_Comm_free(&shm_comm);

    MTest_Finalize(errors);

    if (bases)
        free(bases);

    return MTestReturnValue(errors);
}