source: CPL/oasis3-mct/branches/OASIS3-MCT_5.0_branch/pyoasis/examples/6-apple_and_orange/C/receiver-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: 2.8 KB
Line 
1#include <stdio.h>
2#include <math.h>
3#include "mpi.h"
4#include "oasis_c.h"
5
6int main(int argc, char *argv[])
7{
8  char *comp_name = "receiver";
9  fprintf(stdout,"Component name: %s\n", comp_name);
10  fflush(stdout);
11
12  int comp_id;
13
14  OASIS_CHECK_ERR(oasis_c_init_comp(&comp_id, comp_name, OASIS_COUPLED));
15  fprintf(stdout, "Receiver: Component ID: %d\n", comp_id);
16  fflush(stdout);
17
18  MPI_Comm local_comm;
19  OASIS_CHECK_ERR(oasis_c_get_localcomm(&local_comm));
20  int comm_size;
21  OASIS_CHECK_MPI_ERR(MPI_Comm_size(local_comm, &comm_size));
22  int comm_rank;
23  OASIS_CHECK_MPI_ERR(MPI_Comm_rank(local_comm, &comm_rank));
24  fprintf(stdout,"Receiver: local_comm_rank = %d of %d\n",comm_rank,comm_size);
25  fflush(stdout);
26
27  const int n_points = 16;
28  if ( ( n_points % comm_size ) != 0 ) {
29    OASIS_CHECK_ERR(oasis_c_abort(comp_id, comp_name,
30                                  "Receiver: comm_size has to divide n_points exactly",
31                                  __FILE__, __LINE__));
32  }
33
34  int  local_size = n_points/comm_size;
35  int  offset = comm_rank*local_size;
36
37  int n_segments = 1;
38  int part_id;
39  int i;
40  {
41    int part_params[OASIS_Orange_Params(n_segments)];
42    part_params[OASIS_Strategy] = OASIS_Orange;
43    part_params[OASIS_Segments] = n_segments;
44    for (i = 0; i<n_segments ; i++) {
45      part_params[OASIS_Segments + 2*i + 1] = offset;
46      part_params[OASIS_Segments + 2*i + 2] = local_size;
47    }
48    OASIS_CHECK_ERR(oasis_c_def_partition(&part_id, OASIS_Orange_Params(n_segments),
49                                          part_params, OASIS_No_Gsize,
50                                          OASIS_No_Name));
51    fprintf(stdout, "Receiver rank(%d): part_id: %d\n", comm_rank, part_id);
52    fflush(stdout);
53  }
54
55  char *var_name  = "FRECVATM";
56  fprintf(stdout, "Receiver rank(%d): var_name %s\n", comm_rank, var_name);
57  int bundle_size = 1;
58  int var_id;
59
60  OASIS_CHECK_ERR(oasis_c_def_var(&var_id, var_name, part_id, bundle_size, OASIS_IN, OASIS_REAL));
61  fprintf(stdout, "Receiver rank(%d): var_id %d\n", comm_rank, var_id);
62  fflush(stdout);
63
64  OASIS_CHECK_ERR(oasis_c_enddef());
65
66  float field[local_size];
67  for (i = 0; i<local_size; i++) {
68    field[i] = 0.;
69  }
70  int date = 0;
71
72  int kinfo;
73  OASIS_CHECK_ERR(oasis_c_get(var_id, date, local_size, 1, bundle_size, OASIS_REAL, OASIS_COL_MAJOR, field, &kinfo));
74  fprintf(stdout, "Receiver rank(%d): oasis_c_get returned kinfo = %d\n", comm_rank, kinfo);
75  fflush(stdout);
76
77  OASIS_CHECK_ERR(oasis_c_terminate());
78
79  float epsilon = 1.e-8;
80  float error = 0.;
81
82  for (i = 0; i<local_size; i++) {
83    error += fabs(field[i] - (float) (i+offset));
84  }
85  if (error < epsilon) {
86    fprintf(stdout, "Receiver rank(%d): Data received successfully\n", comm_rank);
87    fflush(stdout);
88  } else {
89    fprintf(stdout, "Receiver rank(%d): Got first four elements\n",comm_rank);
90    for (i = 0; i<4 ; i++ ) {
91      fprintf(stdout, "[rank %d] Element %d contains %f instead of %f\n",
92              comm_rank, i, field[i], (float) i);
93    }
94    fflush(stdout);
95  }
96}
Note: See TracBrowser for help on using the repository browser.