Changeset 27374


Ignore:
Timestamp:
11/10/22 07:10:55 (2 years ago)
Author:
Mathieu Morlighem
Message:

NEW: added non-blocking send/recv

File:
1 edited

Legend:

Unmodified
Added
Removed
  • issm/trunk-jpl/src/c/classes/FemModel.cpp

    r27308 r27374  
    14451445
    14461446        /*Now send and receive vector for vertices on partition edge*/
    1447         #ifdef _HAVE_AD_
    1448         IssmDouble* buffer = xNew<IssmDouble>(this->vertices->Size(),"t"); //only one alloc, "t" is required by adolc
    1449         #else
    1450         IssmDouble* buffer = xNew<IssmDouble>(this->vertices->Size());
    1451         #endif
     1447        IssmDouble       **send_buffers  = xNew<IssmDouble*>(num_procs);
     1448        ISSM_MPI_Request  *send_requests = xNew<ISSM_MPI_Request>(num_procs);
     1449        for (int rank = 0;rank<num_procs;rank++){
     1450                #ifdef _HAVE_AD_
     1451                send_buffers[rank] = xNew<IssmDouble>(this->vertices->Size(),"t"); //only one alloc, "t" is required by adolc
     1452                #else
     1453                send_buffers[rank] = xNew<IssmDouble>(this->vertices->Size());
     1454                #endif
     1455                send_requests[rank] = MPI_REQUEST_NULL;
     1456        }
     1457        IssmDouble* recv_buffer = xNew<IssmDouble>(this->vertices->Size());
    14521458        for(int rank=0;rank<num_procs;rank++){
    14531459                if(this->vertices->common_send[rank]){
     
    14571463                                Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid));
    14581464                                _assert_(!vertex->clone);
    1459                                 buffer[i] = local_vector[vertex->lid];
    1460                         }
    1461                         ISSM_MPI_Send(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm());
     1465            send_buffers[rank][i] = local_vector[vertex->lid];
     1466                        }
     1467         ISSM_MPI_Isend(send_buffers[rank],numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&send_requests[rank]);
    14621468                }
    14631469        }
     
    14651471                if(this->vertices->common_recv[rank]){
    14661472                        int  numids = this->vertices->common_recv[rank];
    1467                         ISSM_MPI_Recv(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
     1473         ISSM_MPI_Recv(recv_buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
    14681474                        for(int i=0;i<numids;i++){
    14691475                                int   master_lid = this->vertices->common_recv_ids[rank][i];
    14701476                                Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid));
    14711477                                _assert_(vertex->clone);
    1472                                 local_vector[vertex->lid] = buffer[i];
    1473                         }
    1474                 }
    1475         }
    1476         xDelete<IssmDouble>(buffer);
     1478            local_vector[vertex->lid] = recv_buffer[i];
     1479                        }
     1480                }
     1481        }
     1482   xDelete<IssmDouble>(recv_buffer);
     1483   for(int rank=0;rank<num_procs;rank++){
     1484      ISSM_MPI_Wait(&send_requests[rank],&status);
     1485      xDelete<IssmDouble>(send_buffers[rank]);
     1486   }
     1487   xDelete<IssmDouble*>(send_buffers);
     1488   xDelete<ISSM_MPI_Request>(send_requests);
    14771489}/*}}}*/
    14781490void FemModel::SyncLocalVectorWithClonesVerticesAdd(IssmDouble* local_vector){/*{{{*/
     
    14841496
    14851497        /*Now send and receive vector for vertices on partition edge*/
    1486         #ifdef _HAVE_AD_
    1487         IssmDouble* buffer = xNew<IssmDouble>(this->vertices->Size(),"t"); //only one alloc, "t" is required by adolc
    1488         #else
    1489         IssmDouble* buffer = xNew<IssmDouble>(this->vertices->Size());
    1490         #endif
     1498        IssmDouble       **send_buffers  = xNew<IssmDouble*>(num_procs);
     1499        ISSM_MPI_Request  *send_requests = xNew<ISSM_MPI_Request>(num_procs);
     1500        for (int rank = 0;rank<num_procs;rank++){
     1501                #ifdef _HAVE_AD_
     1502                send_buffers[rank] = xNew<IssmDouble>(this->vertices->Size(),"t"); //only one alloc, "t" is required by adolc
     1503                #else
     1504                send_buffers[rank] = xNew<IssmDouble>(this->vertices->Size());
     1505                #endif
     1506                send_requests[rank] = MPI_REQUEST_NULL;
     1507        }
     1508        IssmDouble* recv_buffer = xNew<IssmDouble>(this->vertices->Size());
    14911509
    14921510        /*1st: add slaves to master values (reverse of what we usually do)*/
     
    14981516                                Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid));
    14991517                                _assert_(vertex->clone);
    1500                                 buffer[i] = local_vector[vertex->lid];
    1501                         }
    1502                         ISSM_MPI_Send(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm());
     1518                                send_buffers[rank][i] = local_vector[vertex->lid];
     1519                        }
     1520                        ISSM_MPI_Isend(send_buffers[rank],numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&send_requests[rank]);
    15031521                }
    15041522        }
     
    15061524                if(this->vertices->common_send[rank]){
    15071525                        int  numids = this->vertices->common_send[rank];
    1508                         ISSM_MPI_Recv(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
     1526                        ISSM_MPI_Recv(recv_buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
    15091527                        for(int i=0;i<numids;i++){
    15101528                                int   master_lid = this->vertices->common_send_ids[rank][i];
    15111529                                Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid));
    15121530                                _assert_(!vertex->clone);
    1513                                 local_vector[vertex->lid] += buffer[i];
    1514                         }
    1515                 }
    1516         }
     1531                                local_vector[vertex->lid] += recv_buffer[i];
     1532                        }
     1533                }
     1534        }
     1535
     1536        /*Wait until MPI is done*/
     1537        for(int rank=0;rank<num_procs;rank++) ISSM_MPI_Wait(&send_requests[rank],&status);
    15171538
    15181539        /*Now sync masters across partitions*/
     
    15241545                                Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid));
    15251546                                _assert_(!vertex->clone);
    1526                                 buffer[i] = local_vector[vertex->lid];
    1527                         }
    1528                         ISSM_MPI_Send(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm());
     1547                                send_buffers[rank][i] = local_vector[vertex->lid];
     1548                        }
     1549                        ISSM_MPI_Isend(send_buffers[rank],numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&send_requests[rank]);
    15291550                }
    15301551        }
     
    15321553                if(this->vertices->common_recv[rank]){
    15331554                        int  numids = this->vertices->common_recv[rank];
    1534                         ISSM_MPI_Recv(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
     1555                        ISSM_MPI_Recv(recv_buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
     1556
    15351557                        for(int i=0;i<numids;i++){
    15361558                                int   master_lid = this->vertices->common_recv_ids[rank][i];
    15371559                                Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid));
    15381560                                _assert_(vertex->clone);
    1539                                 local_vector[vertex->lid] = buffer[i];
    1540                         }
    1541                 }
    1542         }
    1543         xDelete<IssmDouble>(buffer);
     1561                                local_vector[vertex->lid] = recv_buffer[i];
     1562                        }
     1563                }
     1564        }
     1565        xDelete<IssmDouble>(recv_buffer);
     1566        for(int rank=0;rank<num_procs;rank++){
     1567                ISSM_MPI_Wait(&send_requests[rank],&status);
     1568                xDelete<IssmDouble>(send_buffers[rank]);
     1569        }
     1570        xDelete<IssmDouble*>(send_buffers);
     1571        xDelete<ISSM_MPI_Request>(send_requests);
    15441572}/*}}}*/
    15451573void FemModel::GetLocalVectorWithClonesNodes(IssmDouble** plocal_vector,Vector<IssmDouble> *vector){/*{{{*/
     
    15671595
    15681596        /*Now send and receive vector for nodes on partition edge*/
    1569         #ifdef _HAVE_AD_
    1570         IssmDouble* buffer = xNew<IssmDouble>(this->nodes->Size(),"t"); //only one alloc, "t" is required by adolc
    1571         #else
    1572         IssmDouble* buffer = xNew<IssmDouble>(this->nodes->Size());
    1573         #endif
     1597        /*Now send and receive vector for vertices on partition edge*/
     1598        IssmDouble       **send_buffers  = xNew<IssmDouble*>(num_procs);
     1599        ISSM_MPI_Request  *send_requests = xNew<ISSM_MPI_Request>(num_procs);
     1600        for (int rank = 0;rank<num_procs;rank++){
     1601                #ifdef _HAVE_AD_
     1602                send_buffers[rank] = xNew<IssmDouble>(this->vertices->Size(),"t"); //only one alloc, "t" is required by adolc
     1603                #else
     1604                send_buffers[rank] = xNew<IssmDouble>(this->vertices->Size());
     1605                #endif
     1606                send_requests[rank] = MPI_REQUEST_NULL;
     1607        }
     1608        IssmDouble* recv_buffer = xNew<IssmDouble>(this->vertices->Size());
     1609
    15741610        for(int rank=0;rank<num_procs;rank++){
    15751611                if(this->nodes->common_send[rank]){
     
    15791615                                Node* vertex=xDynamicCast<Node*>(this->nodes->GetObjectByOffset(master_lid));
    15801616                                _assert_(!vertex->clone);
    1581                                 buffer[i] = local_vector[vertex->lid];
    1582                         }
    1583                         ISSM_MPI_Send(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm());
     1617                                send_buffers[rank][i] = local_vector[vertex->lid];
     1618                        }
     1619                        ISSM_MPI_Isend(send_buffers[rank],numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&send_requests[rank]);
    15841620                }
    15851621        }
     
    15871623                if(this->nodes->common_recv[rank]){
    15881624                        int  numids = this->nodes->common_recv[rank];
    1589                         ISSM_MPI_Recv(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
     1625                        ISSM_MPI_Recv(recv_buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
    15901626                        for(int i=0;i<numids;i++){
    15911627                                int   master_lid = this->nodes->common_recv_ids[rank][i];
    15921628                                Node* vertex=xDynamicCast<Node*>(this->nodes->GetObjectByOffset(master_lid));
    15931629                                _assert_(vertex->clone);
    1594                                 local_vector[vertex->lid] = buffer[i];
    1595                         }
    1596                 }
    1597         }
    1598         xDelete<IssmDouble>(buffer);
     1630                                local_vector[vertex->lid] = recv_buffer[i];
     1631                        }
     1632                }
     1633        }
     1634
     1635        xDelete<IssmDouble>(recv_buffer);
     1636        for(int rank=0;rank<num_procs;rank++){
     1637                ISSM_MPI_Wait(&send_requests[rank],&status);
     1638                xDelete<IssmDouble>(send_buffers[rank]);
     1639        }
     1640        xDelete<IssmDouble*>(send_buffers);
     1641        xDelete<ISSM_MPI_Request>(send_requests);
    15991642
    16001643        /*Assign output pointer*/
Note: See TracChangeset for help on using the changeset viewer.