xref: /petsc/src/vec/is/sf/impls/basic/allgatherv/sfallgatherv.c (revision d7cc930e14e615e9907267aaa472dd0ccceeab82)
1 #include <../src/vec/is/sf/impls/basic/allgatherv/sfallgatherv.h>
2 
3 PETSC_INTERN PetscErrorCode PetscSFBcastAndOpBegin_Gatherv(PetscSF,MPI_Datatype,PetscMemType,const void*,PetscMemType,void*,MPI_Op);
4 
5 /* PetscSFGetGraph is non-collective. An implementation should not have collective calls */
6 PETSC_INTERN PetscErrorCode PetscSFGetGraph_Allgatherv(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote)
7 {
8   PetscErrorCode ierr;
9   PetscInt       i,j,k;
10   const PetscInt *range;
11   PetscMPIInt    size;
12 
13   PetscFunctionBegin;
14   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
15   if (nroots)  *nroots  = sf->nroots;
16   if (nleaves) *nleaves = sf->nleaves;
17   if (ilocal)  *ilocal  = NULL; /* Contiguous leaves */
18   if (iremote) {
19     if (!sf->remote && sf->nleaves) { /* The && sf->nleaves makes sfgatherv able to inherit this routine */
20       ierr = PetscLayoutGetRanges(sf->map,&range);CHKERRQ(ierr);
21       ierr = PetscMalloc1(sf->nleaves,&sf->remote);CHKERRQ(ierr);
22       sf->remote_alloc = sf->remote;
23       for (i=0; i<size; i++) {
24         for (j=range[i],k=0; j<range[i+1]; j++,k++) {
25           sf->remote[j].rank  = i;
26           sf->remote[j].index = k;
27         }
28       }
29     }
30     *iremote = sf->remote;
31   }
32   PetscFunctionReturn(0);
33 }
34 
35 PETSC_INTERN PetscErrorCode PetscSFSetUp_Allgatherv(PetscSF sf)
36 {
37   PetscErrorCode     ierr;
38   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv*)sf->data;
39   PetscMPIInt        size;
40   PetscInt           i;
41   const PetscInt     *range;
42 
43   PetscFunctionBegin;
44   ierr = PetscSFSetUp_Allgather(sf);CHKERRQ(ierr);
45   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
46   if (sf->nleaves) { /* This if (sf->nleaves) test makes sfgatherv able to inherit this routine */
47     ierr = PetscMalloc1(size,&dat->recvcounts);CHKERRQ(ierr);
48     ierr = PetscMalloc1(size,&dat->displs);CHKERRQ(ierr);
49     ierr = PetscLayoutGetRanges(sf->map,&range);CHKERRQ(ierr);
50 
51     for (i=0; i<size; i++) {
52       ierr = PetscMPIIntCast(range[i],&dat->displs[i]);CHKERRQ(ierr);
53       ierr = PetscMPIIntCast(range[i+1]-range[i],&dat->recvcounts[i]);CHKERRQ(ierr);
54     }
55   }
56   PetscFunctionReturn(0);
57 }
58 
59 PETSC_INTERN PetscErrorCode PetscSFReset_Allgatherv(PetscSF sf)
60 {
61   PetscErrorCode         ierr;
62   PetscSF_Allgatherv     *dat = (PetscSF_Allgatherv*)sf->data;
63 
64   PetscFunctionBegin;
65   ierr = PetscFree(dat->iranks);CHKERRQ(ierr);
66   ierr = PetscFree(dat->ioffset);CHKERRQ(ierr);
67   ierr = PetscFree(dat->irootloc);CHKERRQ(ierr);
68   ierr = PetscFree(dat->recvcounts);CHKERRQ(ierr);
69   ierr = PetscFree(dat->displs);CHKERRQ(ierr);
70   if (dat->inuse) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Outstanding operation has not been completed");
71   ierr = PetscSFLinkDestroy(sf,&dat->avail);CHKERRQ(ierr);
72   PetscFunctionReturn(0);
73 }
74 
75 PETSC_INTERN PetscErrorCode PetscSFDestroy_Allgatherv(PetscSF sf)
76 {
77   PetscErrorCode ierr;
78 
79   PetscFunctionBegin;
80   ierr = PetscSFReset_Allgatherv(sf);CHKERRQ(ierr);
81   ierr = PetscFree(sf->data);CHKERRQ(ierr);
82   PetscFunctionReturn(0);
83 }
84 
85 static PetscErrorCode PetscSFBcastAndOpBegin_Allgatherv(PetscSF sf,MPI_Datatype unit,PetscMemType rootmtype,const void *rootdata,PetscMemType leafmtype,void *leafdata,MPI_Op op)
86 {
87   PetscErrorCode         ierr;
88   PetscSFLink            link;
89   PetscMPIInt            sendcount;
90   MPI_Comm               comm;
91   void                   *rootbuf = NULL,*leafbuf = NULL;
92   MPI_Request            *req;
93   PetscSF_Allgatherv     *dat = (PetscSF_Allgatherv*)sf->data;
94 
95   PetscFunctionBegin;
96   ierr = PetscSFLinkCreate(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op,PETSCSF_BCAST,&link);CHKERRQ(ierr);
97   ierr = PetscSFLinkPackRootData(sf,link,PETSCSF_REMOTE,rootdata);CHKERRQ(ierr);
98   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
99   ierr = PetscMPIIntCast(sf->nroots,&sendcount);CHKERRQ(ierr);
100   ierr = PetscSFLinkGetMPIBuffersAndRequests(sf,link,PETSCSF_ROOT2LEAF,&rootbuf,&leafbuf,&req,NULL);CHKERRQ(ierr);
101   ierr = MPIU_Iallgatherv(rootbuf,sendcount,unit,leafbuf,dat->recvcounts,dat->displs,unit,comm,req);CHKERRQ(ierr);
102   PetscFunctionReturn(0);
103 }
104 
105 static PetscErrorCode PetscSFReduceBegin_Allgatherv(PetscSF sf,MPI_Datatype unit,PetscMemType leafmtype,const void *leafdata,PetscMemType rootmtype,void *rootdata,MPI_Op op)
106 {
107   PetscErrorCode         ierr;
108   PetscSFLink            link;
109   PetscSF_Allgatherv     *dat = (PetscSF_Allgatherv*)sf->data;
110   PetscInt               rstart;
111   PetscMPIInt            rank,count,recvcount;
112   MPI_Comm               comm;
113   void                   *rootbuf = NULL,*leafbuf = NULL;
114   MPI_Request            *req;
115 
116   PetscFunctionBegin;
117   ierr = PetscSFLinkCreate(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op,PETSCSF_REDUCE,&link);CHKERRQ(ierr);
118   if (op == MPIU_REPLACE) {
119     /* REPLACE is only meaningful when all processes have the same leafdata to reduce. Therefore copying from local leafdata is fine */
120     ierr = PetscLayoutGetRange(sf->map,&rstart,NULL);CHKERRQ(ierr);
121     ierr = (*link->Memcpy)(link,rootmtype,rootdata,leafmtype,(const char*)leafdata+(size_t)rstart*link->unitbytes,(size_t)sf->nroots*link->unitbytes);CHKERRQ(ierr);
122   } else {
123     /* Reduce leafdata, then scatter to rootdata */
124     ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
125     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
126     ierr = PetscSFLinkPackLeafData(sf,link,PETSCSF_REMOTE,leafdata);CHKERRQ(ierr);
127     ierr = PetscSFLinkGetMPIBuffersAndRequests(sf,link,PETSCSF_LEAF2ROOT,&rootbuf,&leafbuf,&req,NULL);CHKERRQ(ierr);
128     ierr = PetscMPIIntCast(dat->rootbuflen[PETSCSF_REMOTE],&recvcount);CHKERRQ(ierr);
129     /* Allocate a separate leaf buffer on rank 0 */
130     if (!rank && !link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi]) {
131       ierr = PetscSFMalloc(sf,link->leafmtype_mpi,sf->leafbuflen[PETSCSF_REMOTE]*link->unitbytes,(void**)&link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi]);CHKERRQ(ierr);
132     }
133     /* In case we already copied leafdata from device to host (i.e., no use_gpu_aware_mpi), we need to adjust leafbuf on rank 0 */
134     if (!rank && link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi] == leafbuf) leafbuf = MPI_IN_PLACE;
135     ierr = PetscMPIIntCast(sf->nleaves*link->bs,&count);CHKERRQ(ierr);
136     ierr = MPI_Reduce(leafbuf,link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi],count,link->basicunit,op,0,comm);CHKERRQ(ierr); /* Must do reduce with MPI builltin datatype basicunit */
137     ierr = MPIU_Iscatterv(link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi],dat->recvcounts,dat->displs,unit,rootbuf,recvcount,unit,0,comm,req);CHKERRQ(ierr);
138   }
139   PetscFunctionReturn(0);
140 }
141 
142 static PetscErrorCode PetscSFBcastToZero_Allgatherv(PetscSF sf,MPI_Datatype unit,PetscMemType rootmtype,const void *rootdata,PetscMemType leafmtype,void *leafdata)
143 {
144   PetscErrorCode         ierr;
145   PetscSFLink            link;
146   PetscMPIInt            rank;
147 
148   PetscFunctionBegin;
149   ierr = PetscSFBcastAndOpBegin_Gatherv(sf,unit,rootmtype,rootdata,leafmtype,leafdata,MPIU_REPLACE);CHKERRQ(ierr);
150   ierr = PetscSFLinkGetInUse(sf,unit,rootdata,leafdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
151   ierr = PetscSFLinkMPIWaitall(sf,link,PETSCSF_ROOT2LEAF);CHKERRQ(ierr);
152   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
153   if (!rank && leafmtype == PETSC_MEMTYPE_DEVICE && !sf->use_gpu_aware_mpi) {
154     ierr = (*link->Memcpy)(link,PETSC_MEMTYPE_DEVICE,leafdata,PETSC_MEMTYPE_HOST,link->leafbuf[PETSC_MEMTYPE_HOST],sf->leafbuflen[PETSCSF_REMOTE]*link->unitbytes);CHKERRQ(ierr);
155   }
156   ierr = PetscSFLinkReclaim(sf,&link);CHKERRQ(ierr);
157   PetscFunctionReturn(0);
158 }
159 
160 /* This routine is very tricky (I believe it is rarely used with this kind of graph so just provide a simple but not-optimal implementation).
161 
162    Suppose we have three ranks. Rank 0 has a root with value 1. Rank 0,1,2 has a leaf with value 2,3,4 respectively. The leaves are connected
163    to the root on rank 0. Suppose op=MPI_SUM and rank 0,1,2 gets root state in their rank order. By definition of this routine, rank 0 sees 1
164    in root, fetches it into its leafupate, then updates root to 1 + 2 = 3; rank 1 sees 3 in root, fetches it into its leafupate, then updates
165    root to 3 + 3 = 6; rank 2 sees 6 in root, fetches it into its leafupdate, then updates root to 6 + 4 = 10.  At the end, leafupdate on rank
166    0,1,2 is 1,3,6 respectively. root is 10.
167 
168    We use a simpler implementation. From the same initial state, we copy leafdata to leafupdate
169              rank-0   rank-1    rank-2
170         Root     1
171         Leaf     2       3         4
172      Leafupdate  2       3         4
173 
174    Do MPI_Exscan on leafupdate,
175              rank-0   rank-1    rank-2
176         Root     1
177         Leaf     2       3         4
178      Leafupdate  2       2         5
179 
180    BcastAndOp from root to leafupdate,
181              rank-0   rank-1    rank-2
182         Root     1
183         Leaf     2       3         4
184      Leafupdate  3       3         6
185 
186    Copy root to leafupdate on rank-0
187              rank-0   rank-1    rank-2
188         Root     1
189         Leaf     2       3         4
190      Leafupdate  1       3         6
191 
192    Reduce from leaf to root,
193              rank-0   rank-1    rank-2
194         Root     10
195         Leaf     2       3         4
196      Leafupdate  1       3         6
197 */
198 PETSC_INTERN PetscErrorCode PetscSFFetchAndOpBegin_Allgatherv(PetscSF sf,MPI_Datatype unit,PetscMemType rootmtype,void *rootdata,PetscMemType leafmtype,const void *leafdata,void *leafupdate,MPI_Op op)
199 {
200   PetscErrorCode         ierr;
201   PetscSFLink            link;
202   MPI_Comm               comm;
203   PetscMPIInt            count;
204 
205   PetscFunctionBegin;
206   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
207   if (!sf->use_gpu_aware_mpi && (rootmtype == PETSC_MEMTYPE_DEVICE || leafmtype == PETSC_MEMTYPE_DEVICE)) SETERRQ(comm,PETSC_ERR_SUP,"Do FetchAndOp on device but not use gpu-aware MPI");
208   /* Copy leafdata to leafupdate */
209   ierr = PetscSFLinkCreate(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op,PETSCSF_FETCH,&link);CHKERRQ(ierr);
210   ierr = PetscSFLinkPackLeafData(sf,link,PETSCSF_REMOTE,leafdata);CHKERRQ(ierr); /* Sync the device */
211   ierr = (*link->Memcpy)(link,leafmtype,leafupdate,leafmtype,leafdata,sf->nleaves*link->unitbytes);CHKERRQ(ierr);
212   ierr = PetscSFLinkGetInUse(sf,unit,rootdata,leafdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr);
213 
214   /* Exscan on leafupdate and then BcastAndOp rootdata to leafupdate */
215   if (op == MPIU_REPLACE) {
216     PetscMPIInt size,rank,prev,next;
217     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
218     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
219     prev = rank ?            rank-1 : MPI_PROC_NULL;
220     next = (rank < size-1) ? rank+1 : MPI_PROC_NULL;
221     ierr = PetscMPIIntCast(sf->nleaves,&count);CHKERRQ(ierr);
222     ierr = MPI_Sendrecv_replace(leafupdate,count,unit,next,link->tag,prev,link->tag,comm,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
223   } else {
224     ierr = PetscMPIIntCast(sf->nleaves*link->bs,&count);CHKERRQ(ierr);
225     ierr = MPI_Exscan(MPI_IN_PLACE,leafupdate,count,link->basicunit,op,comm);CHKERRQ(ierr);
226   }
227   ierr = PetscSFLinkReclaim(sf,&link);CHKERRQ(ierr);
228   ierr = PetscSFBcastAndOpBegin(sf,unit,rootdata,leafupdate,op);CHKERRQ(ierr);
229   ierr = PetscSFBcastAndOpEnd(sf,unit,rootdata,leafupdate,op);CHKERRQ(ierr);
230 
231   /* Bcast roots to rank 0's leafupdate */
232   ierr = PetscSFBcastToZero_Private(sf,unit,rootdata,leafupdate);CHKERRQ(ierr); /* Using this line makes Allgather SFs able to inherit this routine */
233 
234   /* Reduce leafdata to rootdata */
235   ierr = PetscSFReduceBegin(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
236   PetscFunctionReturn(0);
237 }
238 
239 PETSC_INTERN PetscErrorCode PetscSFFetchAndOpEnd_Allgatherv(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
240 {
241   PetscErrorCode         ierr;
242 
243   PetscFunctionBegin;
244   ierr = PetscSFReduceEnd(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
245   PetscFunctionReturn(0);
246 }
247 
248 /* Get root ranks accessing my leaves */
249 PETSC_INTERN PetscErrorCode PetscSFGetRootRanks_Allgatherv(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote)
250 {
251   PetscErrorCode ierr;
252   PetscInt       i,j,k,size;
253   const PetscInt *range;
254 
255   PetscFunctionBegin;
256   /* Lazily construct these large arrays if users really need them for this type of SF. Very likely, they do not */
257   if (sf->nranks && !sf->ranks) { /* On rank!=0, sf->nranks=0. The sf->nranks test makes this routine also works for sfgatherv */
258     size = sf->nranks;
259     ierr = PetscLayoutGetRanges(sf->map,&range);CHKERRQ(ierr);
260     ierr = PetscMalloc4(size,&sf->ranks,size+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);CHKERRQ(ierr);
261     for (i=0; i<size; i++) sf->ranks[i] = i;
262     ierr = PetscArraycpy(sf->roffset,range,size+1);CHKERRQ(ierr);
263     for (i=0; i<sf->nleaves; i++) sf->rmine[i] = i; /*rmine are never NULL even for contiguous leaves */
264     for (i=0; i<size; i++) {
265       for (j=range[i],k=0; j<range[i+1]; j++,k++) sf->rremote[j] = k;
266     }
267   }
268 
269   if (nranks)  *nranks  = sf->nranks;
270   if (ranks)   *ranks   = sf->ranks;
271   if (roffset) *roffset = sf->roffset;
272   if (rmine)   *rmine   = sf->rmine;
273   if (rremote) *rremote = sf->rremote;
274   PetscFunctionReturn(0);
275 }
276 
277 /* Get leaf ranks accessing my roots */
278 PETSC_INTERN PetscErrorCode PetscSFGetLeafRanks_Allgatherv(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc)
279 {
280   PetscErrorCode     ierr;
281   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv*)sf->data;
282   MPI_Comm           comm;
283   PetscMPIInt        size,rank;
284   PetscInt           i,j;
285 
286   PetscFunctionBegin;
287   /* Lazily construct these large arrays if users really need them for this type of SF. Very likely, they do not */
288   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
289   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
290   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
291   if (niranks) *niranks = size;
292 
293   /* PetscSF_Basic has distinguished incoming ranks. Here we do not need that. But we must put self as the first and
294      sort other ranks. See comments in PetscSFSetUp_Basic about MatGetBrowsOfAoCols_MPIAIJ on why.
295    */
296   if (iranks) {
297     if (!dat->iranks) {
298       ierr = PetscMalloc1(size,&dat->iranks);CHKERRQ(ierr);
299       dat->iranks[0] = rank;
300       for (i=0,j=1; i<size; i++) {if (i == rank) continue; dat->iranks[j++] = i;}
301     }
302     *iranks = dat->iranks; /* dat->iranks was init'ed to NULL by PetscNewLog */
303   }
304 
305   if (ioffset) {
306     if (!dat->ioffset) {
307       ierr = PetscMalloc1(size+1,&dat->ioffset);CHKERRQ(ierr);
308       for (i=0; i<=size; i++) dat->ioffset[i] = i*sf->nroots;
309     }
310     *ioffset = dat->ioffset;
311   }
312 
313   if (irootloc) {
314     if (!dat->irootloc) {
315       ierr = PetscMalloc1(sf->nleaves,&dat->irootloc);CHKERRQ(ierr);
316       for (i=0; i<size; i++) {
317         for (j=0; j<sf->nroots; j++) dat->irootloc[i*sf->nroots+j] = j;
318       }
319     }
320     *irootloc = dat->irootloc;
321   }
322   PetscFunctionReturn(0);
323 }
324 
325 PETSC_INTERN PetscErrorCode PetscSFCreateLocalSF_Allgatherv(PetscSF sf,PetscSF *out)
326 {
327   PetscInt       i,nroots,nleaves,rstart,*ilocal;
328   PetscSFNode    *iremote;
329   PetscSF        lsf;
330   PetscErrorCode ierr;
331 
332   PetscFunctionBegin;
333   nleaves = sf->nleaves ? sf->nroots : 0; /* sf->nleaves can be zero with SFGather(v) */
334   nroots  = nleaves;
335   ierr    = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr);
336   ierr    = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr);
337   ierr    = PetscLayoutGetRange(sf->map,&rstart,NULL);CHKERRQ(ierr);
338 
339   for (i=0; i<nleaves; i++) {
340     ilocal[i]        = rstart + i; /* lsf does not change leave indices */
341     iremote[i].rank  = 0;          /* rank in PETSC_COMM_SELF */
342     iremote[i].index = i;          /* root index */
343   }
344 
345   ierr = PetscSFCreate(PETSC_COMM_SELF,&lsf);CHKERRQ(ierr);
346   ierr = PetscSFSetGraph(lsf,nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
347   ierr = PetscSFSetUp(lsf);CHKERRQ(ierr);
348   *out = lsf;
349   PetscFunctionReturn(0);
350 }
351 
352 PETSC_INTERN PetscErrorCode PetscSFCreate_Allgatherv(PetscSF sf)
353 {
354   PetscErrorCode     ierr;
355   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv*)sf->data;
356 
357   PetscFunctionBegin;
358   sf->ops->BcastAndOpEnd   = PetscSFBcastAndOpEnd_Basic;
359   sf->ops->ReduceEnd       = PetscSFReduceEnd_Basic;
360 
361   sf->ops->SetUp           = PetscSFSetUp_Allgatherv;
362   sf->ops->Reset           = PetscSFReset_Allgatherv;
363   sf->ops->Destroy         = PetscSFDestroy_Allgatherv;
364   sf->ops->GetRootRanks    = PetscSFGetRootRanks_Allgatherv;
365   sf->ops->GetLeafRanks    = PetscSFGetLeafRanks_Allgatherv;
366   sf->ops->GetGraph        = PetscSFGetGraph_Allgatherv;
367   sf->ops->BcastAndOpBegin = PetscSFBcastAndOpBegin_Allgatherv;
368   sf->ops->ReduceBegin     = PetscSFReduceBegin_Allgatherv;
369   sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Allgatherv;
370   sf->ops->FetchAndOpEnd   = PetscSFFetchAndOpEnd_Allgatherv;
371   sf->ops->CreateLocalSF   = PetscSFCreateLocalSF_Allgatherv;
372   sf->ops->BcastToZero     = PetscSFBcastToZero_Allgatherv;
373 
374   ierr = PetscNewLog(sf,&dat);CHKERRQ(ierr);
375   sf->data = (void*)dat;
376   PetscFunctionReturn(0);
377 }
378