Changeset 23642 for issm/trunk-jpl/src/c/classes/FemModel.cpp
- Timestamp:
- 01/17/19 12:57:57 (6 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk-jpl/src/c/classes/FemModel.cpp
r23638 r23642 1422 1422 *plocal_vector = local_vector; 1423 1423 }/*}}}*/ 1424 void FemModel::GetLocalVectorWithClonesNodes(IssmDouble** plocal_vector,Vector<IssmDouble> *vector){/*{{{*/ 1425 1426 /*recover my_rank:*/ 1427 ISSM_MPI_Status status; 1428 int my_rank = IssmComm::GetRank(); 1429 int num_procs = IssmComm::GetSize(); 1430 1431 /*retrieve vertex info*/ 1432 int localsize = this->nodes->NumberOfNodesLocalAll(); 1433 int localsize_masters = this->nodes->NumberOfNodesLocal(); 1434 1435 /*Get local vector of vector*/ 1436 int *indices_vector_masters = NULL; 1437 IssmDouble *local_vector_masters = NULL; 1438 vector->GetLocalVector(&local_vector_masters,&indices_vector_masters); 1439 _assert_(localsize_masters==indices_vector_masters[localsize_masters-1] - indices_vector_masters[0]+1); 1440 xDelete<int>(indices_vector_masters); 1441 1442 /*Now, extend vectors to account for clones (make vectors longer, for clones at the end)*/ 1443 IssmDouble *local_vector = xNew<IssmDouble>(localsize); 1444 xMemCpy<IssmDouble>(local_vector,local_vector_masters,localsize_masters); 1445 xDelete<IssmDouble>(local_vector_masters); 1446 1447 /*Now send and receive vector for nodes on partition edge*/ 1448 #ifdef _HAVE_AD_ 1449 IssmDouble* buffer = xNew<IssmDouble>(this->nodes->Size(),"t"); //only one alloc, "t" is required by adolc 1450 #else 1451 IssmDouble* buffer = xNew<IssmDouble>(this->nodes->Size()); 1452 #endif 1453 for(int rank=0;rank<num_procs;rank++){ 1454 if(this->nodes->common_send[rank]){ 1455 int numids = this->nodes->common_send[rank]; 1456 for(int i=0;i<numids;i++){ 1457 int master_lid = this->nodes->common_send_ids[rank][i]; 1458 Node* vertex=xDynamicCast<Node*>(this->nodes->GetObjectByOffset(master_lid)); 1459 _assert_(!vertex->clone); 1460 buffer[i] = local_vector[vertex->lid]; 1461 } 1462 ISSM_MPI_Send(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm()); 1463 } 1464 } 1465 for(int rank=0;rank<num_procs;rank++){ 1466 if(this->nodes->common_recv[rank]){ 1467 int numids = this->nodes->common_recv[rank]; 1468 ISSM_MPI_Recv(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status); 1469 for(int i=0;i<numids;i++){ 1470 int master_lid = this->nodes->common_recv_ids[rank][i]; 1471 Node* vertex=xDynamicCast<Node*>(this->nodes->GetObjectByOffset(master_lid)); 1472 _assert_(vertex->clone); 1473 local_vector[vertex->lid] = buffer[i]; 1474 } 1475 } 1476 } 1477 xDelete<IssmDouble>(buffer); 1478 1479 /*Assign output pointer*/ 1480 *plocal_vector = local_vector; 1481 }/*}}}*/ 1424 1482 void FemModel::GroundedAreax(IssmDouble* pV, bool scaled){/*{{{*/ 1425 1483
Note:
See TracChangeset
for help on using the changeset viewer.