source: CPL/oasis3-mct/branches/OASIS3-MCT_5.0_branch/pyoasis/examples/4-orange/C/sender-orange.c @ 6331

Last change on this file since 6331 was 6331, checked in by aclsce, 17 months ago

Moved oasis-mct_5.0 in oasis3-mct/branches directory.

File size: 3.0 KB
Line 
1#include <stdio.h>
2#include "mpi.h"
3#include "oasis_c.h"
4
5int main(int argc, char *argv[])
6{
7  char *comp_name = "sender-orange";
8  fprintf(stdout,"Component name: %s\n", comp_name);
9  fflush(stdout);
10
11  int comp_id;
12
13  OASIS_CHECK_ERR(oasis_c_init_comp(&comp_id, comp_name, OASIS_COUPLED));
14  fprintf(stdout, "Sender: Component ID: %d\n", comp_id);
15  fflush(stdout);
16
17  MPI_Comm local_comm;
18  OASIS_CHECK_ERR(oasis_c_get_localcomm(&local_comm));
19  int local_comm_size;
20  OASIS_CHECK_MPI_ERR(MPI_Comm_size(local_comm, &local_comm_size));
21  int local_comm_rank;
22  OASIS_CHECK_MPI_ERR(MPI_Comm_rank(local_comm, &local_comm_rank));
23  fprintf(stdout,"Sender: local_comm_rank = %d of %d\n",local_comm_rank,local_comm_size);
24  fflush(stdout);
25
26  int icpl = 1;
27  if (local_comm_size > 3) {
28    if (local_comm_rank >= local_comm_size - 2) icpl = MPI_UNDEFINED;
29  }
30
31  MPI_Comm coupl_comm;
32  OASIS_CHECK_ERR(oasis_c_create_couplcomm(icpl, local_comm, &coupl_comm));
33
34  int local_size;
35  int offset;
36  int bundle_size;
37  int var_id;
38  int i;
39 
40  if ( icpl == 1 ) {
41    int comm_size;
42    OASIS_CHECK_MPI_ERR(MPI_Comm_size(coupl_comm, &comm_size));
43    int comm_rank;
44    OASIS_CHECK_MPI_ERR(MPI_Comm_rank(coupl_comm, &comm_rank));
45    fprintf(stdout,"Sender: comm_rank = %d of %d\n",comm_rank,comm_size);
46    fflush(stdout);
47
48    const int n_points = 16;
49    if ( ( n_points % comm_size ) != 0 ) {
50      OASIS_CHECK_ERR(oasis_c_abort(comp_id, comp_name,
51                                    "Sender: comm_size has to divide n_points exactly",
52                                    __FILE__, __LINE__));
53    }
54    local_size = n_points/comm_size;
55    offset = comm_rank*local_size;
56    int n_segments = 1;
57    int part_id;
58    {
59      int part_params[OASIS_Orange_Params(n_segments)];
60      part_params[OASIS_Strategy] = OASIS_Orange;
61      part_params[OASIS_Segments] = n_segments;
62      for (i = 0; i<n_segments ; i++) {
63        part_params[OASIS_Segments + 2*i + 1] = offset;
64        part_params[OASIS_Segments + 2*i + 2] = local_size;
65      }
66      OASIS_CHECK_ERR(oasis_c_def_partition(&part_id, OASIS_Orange_Params(n_segments),
67                                            part_params, OASIS_No_Gsize,
68                                            OASIS_No_Name));
69      fprintf(stdout, "Sender rank(%d): part_id: %d\n", local_comm_rank, part_id);
70      fflush(stdout);
71    }
72
73    char *var_name  = "FSENDOCN";
74    fprintf(stdout, "Sender rank(%d): var_name %s\n", local_comm_rank, var_name);
75    bundle_size = 1;
76
77    OASIS_CHECK_ERR(oasis_c_def_var(&var_id, var_name, part_id, bundle_size, OASIS_OUT, OASIS_REAL));
78    fprintf(stdout, "Sender rank(%d): var_id %d\n", local_comm_rank, var_id);
79    fflush(stdout);
80  }
81
82  OASIS_CHECK_ERR(oasis_c_enddef());
83
84  if ( icpl == 1 ) {
85    float field[local_size];
86    for (i = 0; i<local_size; i++) {
87      field[i] = i + offset;
88    }
89
90    int date = 0;
91
92    int kinfo;
93    OASIS_CHECK_ERR(oasis_c_put(var_id, date, local_size, 1, bundle_size, OASIS_REAL, OASIS_COL_MAJOR, field, OASIS_No_Restart, &kinfo));
94    fprintf(stdout, "Sender rank(%d): oasis_c_put returned kinfo = %d\n", local_comm_rank, kinfo);
95    fflush(stdout);
96
97  }
98
99  OASIS_CHECK_ERR(oasis_c_terminate());
100
101}
Note: See TracBrowser for help on using the repository browser.