Changeset 24047
- Timestamp:
- 06/24/19 13:37:54 (6 years ago)
- Location:
- issm/trunk-jpl/src/c
- Files:
-
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk-jpl/src/c/classes/FemModel.cpp
r24035 r24047 1349 1349 IssmDouble* buffer = xNew<IssmDouble>(this->vertices->Size()); 1350 1350 #endif 1351 for(int rank=0;rank<num_procs;rank++){ 1352 if(this->vertices->common_send[rank]){ 1353 int numids = this->vertices->common_send[rank]; 1354 for(int i=0;i<numids;i++){ 1355 int master_lid = this->vertices->common_send_ids[rank][i]; 1356 Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid)); 1357 _assert_(!vertex->clone); 1358 buffer[i] = local_vector[vertex->lid]; 1359 } 1360 ISSM_MPI_Send(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm()); 1361 } 1362 } 1363 for(int rank=0;rank<num_procs;rank++){ 1364 if(this->vertices->common_recv[rank]){ 1365 int numids = this->vertices->common_recv[rank]; 1366 ISSM_MPI_Recv(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status); 1367 for(int i=0;i<numids;i++){ 1368 int master_lid = this->vertices->common_recv_ids[rank][i]; 1369 Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid)); 1370 _assert_(vertex->clone); 1371 local_vector[vertex->lid] = buffer[i]; 1372 } 1373 } 1374 } 1375 xDelete<IssmDouble>(buffer); 1376 }/*}}}*/ 1377 void FemModel::SyncLocalVectorWithClonesVerticesAdd(IssmDouble* local_vector){/*{{{*/ 1378 1379 /*recover my_rank:*/ 1380 ISSM_MPI_Status status; 1381 int my_rank = IssmComm::GetRank(); 1382 int num_procs = IssmComm::GetSize(); 1383 1384 /*Now send and receive vector for vertices on partition edge*/ 1385 #ifdef _HAVE_AD_ 1386 IssmDouble* buffer = xNew<IssmDouble>(this->vertices->Size(),"t"); //only one alloc, "t" is required by adolc 1387 #else 1388 IssmDouble* buffer = xNew<IssmDouble>(this->vertices->Size()); 1389 #endif 1390 1391 /*1st: add slaves to master values (reverse of what we usually do)*/ 1392 for(int rank=0;rank<num_procs;rank++){ 1393 if(this->vertices->common_recv[rank]){ 1394 int numids = this->vertices->common_recv[rank]; 1395 for(int i=0;i<numids;i++){ 1396 int master_lid = this->vertices->common_recv_ids[rank][i]; 1397 Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid)); 1398 _assert_(vertex->clone); 1399 buffer[i] = local_vector[vertex->lid]; 1400 } 1401 ISSM_MPI_Send(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm()); 1402 } 1403 } 1404 for(int rank=0;rank<num_procs;rank++){ 1405 if(this->vertices->common_send[rank]){ 1406 int numids = this->vertices->common_send[rank]; 1407 ISSM_MPI_Recv(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status); 1408 for(int i=0;i<numids;i++){ 1409 int master_lid = this->vertices->common_send_ids[rank][i]; 1410 Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid)); 1411 _assert_(!vertex->clone); 1412 local_vector[vertex->lid] += buffer[i]; 1413 } 1414 } 1415 } 1416 1417 /*Now sync masters across partitions*/ 1351 1418 for(int rank=0;rank<num_procs;rank++){ 1352 1419 if(this->vertices->common_send[rank]){ -
issm/trunk-jpl/src/c/classes/FemModel.h
r24035 r24047 99 99 void GetLocalVectorWithClonesVertices(IssmDouble** plocal_vector,Vector<IssmDouble> *vector); 100 100 void SyncLocalVectorWithClonesVertices(IssmDouble* local_vector); 101 void SyncLocalVectorWithClonesVerticesAdd(IssmDouble* local_vector); 101 102 void GetLocalVectorWithClonesNodes(IssmDouble** plocal_vector,Vector<IssmDouble> *vector); 102 103 void GroundedAreax(IssmDouble* pV, bool scaled); -
issm/trunk-jpl/src/c/modules/KillIcebergsx/KillIcebergsx.cpp
r24041 r24047 42 42 43 43 /*Get local mask from parallel vector*/ 44 femmodel->SyncLocalVectorWithClonesVertices (local_mask);44 femmodel->SyncLocalVectorWithClonesVerticesAdd(local_mask); 45 45 46 46 /*Local iterations on partition*/
Note:
See TracChangeset
for help on using the changeset viewer.