source: CPL/oasis3-mct/branches/OASIS3-MCT_5.0_branch/pyoasis/examples/11-test-interpolation/C/sender-apple.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 <math.h>
3#include "mpi.h"
4#include "netcdf.h"
5#include "oasis_c.h"
6
7int main(int argc, char *argv[])
8{
9  char *comp_name = "sender-apple";
10  int comp_id;
11  OASIS_CHECK_ERR(oasis_c_init_comp(&comp_id, comp_name, OASIS_COUPLED));
12
13  MPI_Comm localcomm;
14  OASIS_CHECK_ERR(oasis_c_get_localcomm(&localcomm));
15  int comm_size;
16  OASIS_CHECK_MPI_ERR(MPI_Comm_size(localcomm, &comm_size));
17  int comm_rank;
18  OASIS_CHECK_MPI_ERR(MPI_Comm_rank(localcomm, &comm_rank));
19  if (comm_rank == 0) {
20    fprintf(stdout,"%s: Component ID: %d\n", comp_name, comp_id);
21  }
22
23  int nx_global = 362;
24  int ny_global = 294;
25  int n_points = nx_global*ny_global;
26
27  int ncid, varid;
28  if (nc_open("grids.nc", 0, &ncid) != NC_NOERR)
29    OASIS_CHECK_ERR(oasis_c_abort(comp_id, comp_name, "Sender: Error in opening grids.nc",
30                                    __FILE__, __LINE__));
31
32  double lon[n_points];
33  double lat[n_points];
34  if (nc_inq_varid(ncid, "nogt.lon", &varid) != NC_NOERR)
35    OASIS_CHECK_ERR(oasis_c_abort(comp_id, comp_name, "Sender: Error in getting nogt.lon id",__FILE__, __LINE__));
36  if (nc_get_var(ncid, varid, lon) != NC_NOERR)
37    OASIS_CHECK_ERR(oasis_c_abort(comp_id, comp_name, "Sender: Error in getting nogt.lon values",__FILE__, __LINE__));
38  if (nc_inq_varid(ncid, "nogt.lat", &varid) != NC_NOERR)
39    OASIS_CHECK_ERR(oasis_c_abort(comp_id, comp_name, "Sender: Error in getting nogt.lat id",__FILE__, __LINE__));
40  if (nc_get_var(ncid, varid, lat) != NC_NOERR)
41    OASIS_CHECK_ERR(oasis_c_abort(comp_id, comp_name, "Sender: Error in getting nogt.lat values",__FILE__, __LINE__));
42  if (nc_close(ncid))
43    OASIS_CHECK_ERR(oasis_c_abort(comp_id, comp_name, "Sender: Error in closing grids.nc",__FILE__, __LINE__));
44
45  int local_size = (int) n_points / comm_size ;
46  int offset = comm_rank * local_size;
47  if ( comm_rank == comm_size - 1) local_size = n_points - offset;
48
49  int part_params[OASIS_Apple_Params];
50  part_params[OASIS_Strategy] = OASIS_Apple;
51  part_params[OASIS_Offset] = offset;
52  part_params[OASIS_Length] = local_size;
53  int part_id;
54  OASIS_CHECK_ERR(oasis_c_def_partition(&part_id, OASIS_Apple_Params,
55                                        part_params, OASIS_No_Gsize,
56                                        OASIS_No_Name));
57
58  char *var_name  = "FSENDANA";
59  int bundle_size = 2;
60  int var_id;
61  OASIS_CHECK_ERR(oasis_c_def_var(&var_id, var_name, part_id, bundle_size, OASIS_OUT, OASIS_DOUBLE));
62
63  OASIS_CHECK_ERR(oasis_c_enddef());
64
65  double dp_conv = atan(1.)/45.0;
66  double bundle[local_size][bundle_size];
67  int i;
68  for (i = 0; i<local_size; i++) {
69    bundle[i][0] = 2.0 + pow(sin(2.0 * lat[offset+i]*dp_conv),4.0) *
70      cos(4.0 * lon[offset+i]*dp_conv);
71    bundle[i][1] = 2.0 - cos(atan(1.)*4.*
72                             (acos(cos(lon[offset+i]*dp_conv)*cos(lat[offset+i]*dp_conv))/
73                              (1.2*atan(1.)*4)));
74  }
75  int date = 0;
76
77  int kinfo;
78  OASIS_CHECK_ERR(oasis_c_put(var_id, date, local_size, 1, bundle_size, OASIS_DOUBLE, OASIS_ROW_MAJOR, bundle, OASIS_No_Restart, &kinfo));
79  fprintf(stdout, "Sender rank(%d): oasis_c_put at time %d returned kinfo = %d\n", comm_rank, date, kinfo);
80  fflush(stdout);
81
82  OASIS_CHECK_ERR(oasis_c_terminate());
83
84}
Note: See TracBrowser for help on using the repository browser.