source: CPL/oasis3-mct/branches/OASIS3-MCT_5.0_branch/pyoasis/examples/3-box/C/sender-box.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.5 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-box";
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 localcomm;
18  OASIS_CHECK_ERR(oasis_c_get_localcomm(&localcomm));
19  int comm_size;
20  OASIS_CHECK_MPI_ERR(MPI_Comm_size(localcomm, &comm_size));
21  int comm_rank;
22  OASIS_CHECK_MPI_ERR(MPI_Comm_rank(localcomm, &comm_rank));
23  fprintf(stdout,"Sender: rank = %d of %d\n",comm_rank,comm_size);
24  fflush(stdout);
25
26  if ( comm_size != 4) {
27    OASIS_CHECK_ERR(oasis_c_abort(comp_id, comp_name,
28                                  "Sender: comm_size has to be 4 for this example",
29                                  __FILE__, __LINE__));
30  }
31
32  const int n_points = 16;
33  int local_size = n_points/comm_size;
34  int offsets[4];
35  offsets[0] = 0;
36  offsets[1] = 2;
37  offsets[2] = 8;
38  offsets[3] = 10;
39
40  int part_params[OASIS_Box_Params];
41  part_params[OASIS_Strategy] = OASIS_Box;
42  part_params[OASIS_Offset]   = offsets[comm_rank];
43  part_params[OASIS_SizeX]    = 2;
44  part_params[OASIS_SizeY]    = 2;
45  part_params[OASIS_LdX]      = 4;
46  int part_id;
47
48  OASIS_CHECK_ERR(oasis_c_def_partition(&part_id, OASIS_Box_Params,
49                                        part_params, OASIS_No_Gsize,
50                                        OASIS_No_Name));
51  fprintf(stdout, "Sender rank(%d): part_id: %d\n", comm_rank, part_id);
52  fflush(stdout);
53
54  char *var_name  = "FSENDOCN";
55  fprintf(stdout, "Sender rank(%d): var_name %s\n", comm_rank, var_name);
56  int bundle_size = 1;
57  int var_id;
58
59  OASIS_CHECK_ERR(oasis_c_def_var(&var_id, var_name, part_id, bundle_size, OASIS_OUT, OASIS_REAL));
60  fprintf(stdout, "Sender rank(%d): var_id %d\n", comm_rank, var_id);
61  fflush(stdout);
62
63  char *var_name2  = "NOTANAME";
64  fprintf(stdout, "Sender rank(%d): var_name %s\n", comm_rank, var_name2);
65  int var_id2;
66  OASIS_CHECK_ERR(oasis_c_def_var(&var_id2, var_name2, part_id, bundle_size, OASIS_OUT, OASIS_REAL));
67  if (var_id2 < 0) {
68    fprintf(stdout, "Sender rank(%d): var_id %d is not active\n", comm_rank, var_id2);
69  } else {
70    fprintf(stdout, "Sender rank(%d): var_id %d\n", comm_rank, var_id2);
71  }   
72  fflush(stdout);
73
74  OASIS_CHECK_ERR(oasis_c_enddef());
75
76  MPI_Comm intra_comm;
77  OASIS_CHECK_ERR(oasis_c_get_intracomm(&intra_comm, "receiver"));
78
79  int intra_size;
80  OASIS_CHECK_MPI_ERR(MPI_Comm_size(intra_comm, &intra_size));
81  int intra_rank;
82  OASIS_CHECK_MPI_ERR(MPI_Comm_rank(intra_comm, &intra_rank));
83  fprintf(stdout,"Sender intra_comm: rank = %d of %d\n",intra_rank,intra_size);
84  fflush(stdout);
85
86  MPI_Comm inter_comm;
87  OASIS_CHECK_ERR(oasis_c_get_intercomm(&inter_comm, "receiver"));
88
89  int inter_size;
90  OASIS_CHECK_MPI_ERR(MPI_Comm_size(inter_comm, &inter_size));
91  int inter_rank;
92  OASIS_CHECK_MPI_ERR(MPI_Comm_rank(inter_comm, &inter_rank));
93  int inter_rsize;
94  OASIS_CHECK_MPI_ERR(MPI_Comm_remote_size(inter_comm, &inter_rsize));
95  fprintf(stdout,"Sender inter_comm: rank = %d of %d Remote size = %d\n",
96          inter_rank,inter_size,inter_rsize);
97  fflush(stdout);
98
99  float field[local_size];
100  field[0] = offsets[comm_rank]+0;
101  field[1] = offsets[comm_rank]+1;
102  field[2] = offsets[comm_rank]+4;
103  field[3] = offsets[comm_rank]+5;
104
105  int date = 0;
106
107  int kinfo;
108  OASIS_CHECK_ERR(oasis_c_put(var_id, date, local_size, 1, bundle_size, OASIS_REAL, OASIS_COL_MAJOR, field, OASIS_No_Restart, &kinfo));
109  fprintf(stdout, "Sender rank(%d): oasis_c_put returned kinfo = %d\n", comm_rank, kinfo);
110  fflush(stdout);
111
112  OASIS_CHECK_ERR(oasis_c_terminate());
113
114}
Note: See TracBrowser for help on using the repository browser.