source: CPL/oasis3-mct/branches/OASIS3-MCT_5.0_branch/lib/mct/mpi-serial/comm.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: 4.0 KB
Line 
1
2#include "mpiP.h"
3
4
5
6/*
7 * Communicators
8 *
9 */
10
11
12
13MPI_Comm mpi_comm_new(void)
14{
15  MPI_Comm chandle;
16  Comm *cptr;
17  static int num=0;
18
19  mpi_alloc_handle(&chandle,(void **) &cptr);
20
21  cptr->sendlist=AP_list_new();
22  cptr->recvlist=AP_list_new();
23
24  cptr->num=num++;
25
26  return(chandle);
27}
28
29
30/*********/
31
32
33FC_FUNC( mpi_comm_free , MPI_COMM_FREE )(int *comm, int *ierror)
34{
35  *ierror=MPI_Comm_free(comm);
36}
37
38
39/*
40 * MPI_Comm_free()
41 *
42 * Note: will NOT free any pending MPI_Request handles
43 * that are allocated...  correct user code should have
44 * already done a Wait or Test to free them.
45 *
46 */
47
48
49int MPI_Comm_free(MPI_Comm *comm)
50{
51  pList sendlist, recvlist;
52  int size;
53  Comm *mycomm;
54
55  mycomm=mpi_handle_to_ptr(*comm);   /* (Comm *)(*comm) */
56
57  sendlist=mycomm->sendlist;
58  recvlist=mycomm->recvlist;
59
60  size=AP_list_size(sendlist);
61  if (size!=0)
62    fprintf(stderr,"MPI_Comm_free: warning: %d pending send reqs\n",
63            size);
64  AP_list_free(sendlist);
65
66
67  size=AP_list_size(recvlist);
68  if (size!=0)
69    fprintf(stderr,"MPI_Comm_free: warning: %d pending receive reqs\n",
70            size);
71  AP_list_free(recvlist);
72
73  mpi_free_handle(*comm);            /* free(mycomm); */
74  *comm=MPI_COMM_NULL;
75
76  return(MPI_SUCCESS);
77}
78
79
80/*********/
81
82
83
84FC_FUNC( mpi_comm_size , MPI_COMM_SIZE )(int *comm, int *size, int *ierror)
85{
86  *ierror=MPI_Comm_size(*comm, size);
87}
88
89
90
91int MPI_Comm_size(MPI_Comm comm, int *size)
92{
93  *size=1;
94
95  return(MPI_SUCCESS);
96}
97
98
99/*********/
100
101
102FC_FUNC( mpi_comm_rank , MPI_COMM_RANK )(int *comm, int *rank, int *ierror)
103{
104  *ierror=MPI_Comm_rank( *comm, rank);
105}
106
107
108int MPI_Comm_rank(MPI_Comm comm, int *rank)
109{
110  *rank=0;
111
112  return(MPI_SUCCESS);
113}
114
115
116
117/*********/
118
119
120FC_FUNC( mpi_comm_dup , MPI_COMM_DUP )(int *comm, int *newcomm, int *ierror)
121{
122
123  *ierror=MPI_Comm_dup( *comm, newcomm);
124
125}
126
127
128int MPI_Comm_dup(MPI_Comm comm, MPI_Comm *newcomm)
129{
130  *newcomm= mpi_comm_new();
131
132#ifdef INFO
133  fflush(stdout);
134  fprintf(stderr,"MPI_Comm_dup: new comm handle=%d\n",*newcomm);
135#endif
136
137  return(MPI_SUCCESS);
138}
139
140
141/*********/
142
143
144int FC_FUNC( mpi_comm_create, MPI_COMM_CREATE)
145     (int *comm, int *group, int *newcomm, int *ierror)
146{
147  *ierror=MPI_Comm_create(*comm,*group,newcomm);
148}
149
150
151
152int MPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm *newcomm)
153{
154  if (group==MPI_GROUP_NULL || group==MPI_GROUP_EMPTY)
155    *newcomm= MPI_COMM_NULL;
156  else
157    *newcomm=mpi_comm_new();
158
159  return(MPI_SUCCESS);
160}
161
162
163
164/*********/
165
166
167FC_FUNC( mpi_comm_split, MPI_COMM_SPLIT )
168     (int *comm, int *color, int *key, int *newcomm, int *ierror)
169{
170  *ierror=MPI_Comm_split(*comm,*color,*key,newcomm);
171
172}
173
174
175
176int MPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm *newcomm)
177{
178  if (color==MPI_UNDEFINED)
179    *newcomm=MPI_COMM_NULL;
180  else
181    *newcomm= mpi_comm_new();
182
183  return(MPI_SUCCESS);
184}
185
186
187/*********/
188
189
190FC_FUNC( mpi_comm_group, MPI_COMM_GROUP )
191     (int *comm, int *group, int *ierror)
192{
193  *ierror= MPI_Comm_group(*comm, group);
194}
195
196
197
198int MPI_Comm_group(MPI_Comm comm, MPI_Group *group)
199{
200  if (comm==MPI_COMM_NULL)
201    *group= MPI_GROUP_NULL;
202  else
203    *group= MPI_GROUP_ONE;
204
205  return(MPI_SUCCESS);
206}
207
208/* Intercomm_create
209 *
210 */
211
212FC_FUNC(mpi_intercomm_create, MPI_INTERCOMM_CREATE)(
213                          int * local_comm, int * local_leader,
214                          int * peer_comm,  int * remote_leader,
215                          int * tag, int * newintercomm, int* ierr)
216{
217  *ierr = MPI_Intercomm_create(*local_comm, *local_leader, *peer_comm,
218                               *remote_leader, *tag, newintercomm);
219}
220
221int MPI_Intercomm_create(MPI_Comm local_comm, int local_leader,
222                          MPI_Comm peer_comm, int remote_leader,
223                          int tag, MPI_Comm *newintercomm)
224{
225  if (local_comm==MPI_COMM_NULL && peer_comm==MPI_COMM_NULL)
226    *newintercomm = MPI_COMM_NULL;
227  else
228    MPI_Comm_dup(MPI_COMM_WORLD, newintercomm);
229
230  return MPI_SUCCESS;
231}
232
233
234/*********/
235
236
237MPI_Comm MPI_Comm_f2c(MPI_Fint comm)
238{
239  /* Comm is an integer handle used both by C and Fortran */
240  return(comm);
241}
242
243
244MPI_Fint MPI_Comm_c2f(MPI_Comm comm)
245{
246  return(comm);
247}
Note: See TracBrowser for help on using the repository browser.