Changeset 23926
- Timestamp:
- 05/20/19 16:33:50 (6 years ago)
- Location:
- issm/trunk-jpl/src/c/classes
- Files:
-
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk-jpl/src/c/classes/FemModel.cpp
r23905 r23926 1308 1308 void FemModel::GetLocalVectorWithClonesGset(IssmDouble** plocal_ug,Vector<IssmDouble> *ug){/*{{{*/ 1309 1309 1310 /*recover my_rank:*/ 1311 ISSM_MPI_Status status; 1312 int my_rank = IssmComm::GetRank(); 1313 int num_procs = IssmComm::GetSize(); 1314 1315 /*retrieve node info*/ 1316 int glocalsize = this->nodes->NumberOfDofsLocalAll(GsetEnum); 1317 int glocalsize_masters = this->nodes->NumberOfDofsLocal(GsetEnum); 1318 int maxdofspernode = this->nodes->MaxNumDofs(GsetEnum); 1319 1320 /*Get local vector of ug*/ 1321 int *indices_ug_masters = NULL; 1322 IssmDouble *local_ug_masters = NULL; 1323 ug->GetLocalVector(&local_ug_masters,&indices_ug_masters); 1324 _assert_(glocalsize_masters==indices_ug_masters[glocalsize_masters-1] - indices_ug_masters[0]+1); 1325 xDelete<int>(indices_ug_masters); 1326 1327 /*Now, extend vectors to account for clones (make vectors longer, for clones at the end)*/ 1328 IssmDouble *local_ug = xNew<IssmDouble>(glocalsize); 1329 xMemCpy<IssmDouble>(local_ug,local_ug_masters,glocalsize_masters); 1330 xDelete<IssmDouble>(local_ug_masters); 1331 1332 /*Now send and receive ug for nodes on partition edge*/ 1333 #ifdef _HAVE_AD_ 1334 IssmDouble* buffer = xNew<IssmDouble>(this->nodes->Size()*maxdofspernode,"t"); //only one alloc, "t" is required by adolc 1335 #else 1336 IssmDouble* buffer = xNew<IssmDouble>(this->nodes->Size()*maxdofspernode); 1337 #endif 1338 for(int rank=0;rank<num_procs;rank++){ 1339 if(this->nodes->common_send[rank]){ 1340 int numids = this->nodes->common_send[rank]; 1341 for(int i=0;i<numids;i++){ 1342 int master_lid = this->nodes->common_send_ids[rank][i]; 1343 Node* node=xDynamicCast<Node*>(this->nodes->GetObjectByOffset(master_lid)); 1344 _assert_(!node->IsClone()); 1345 for(int j=0;j<node->gsize;j++) buffer[i*maxdofspernode+j]=local_ug[node->gdoflist_local[j]]; 1346 } 1347 ISSM_MPI_Send(buffer,numids*maxdofspernode,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm()); 1348 } 1349 } 1350 for(int rank=0;rank<num_procs;rank++){ 1351 if(this->nodes->common_recv[rank]){ 1352 int numids = this->nodes->common_recv[rank]; 1353 ISSM_MPI_Recv(buffer,numids*maxdofspernode,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status); 1354 for(int i=0;i<numids;i++){ 1355 int master_lid = this->nodes->common_recv_ids[rank][i]; 1356 Node* node=xDynamicCast<Node*>(this->nodes->GetObjectByOffset(master_lid)); 1357 for(int j=0;j<node->gsize;j++) local_ug[node->gdoflist_local[j]] = buffer[i*maxdofspernode+j]; 1358 } 1359 } 1360 } 1361 xDelete<IssmDouble>(buffer); 1362 1363 /*Assign output pointer*/ 1364 *plocal_ug = local_ug; 1310 this->nodes->GetLocalVectorWithClonesGset(plocal_ug,ug); 1311 1365 1312 }/*}}}*/ 1366 1313 void FemModel::GetLocalVectorWithClonesVertices(IssmDouble** plocal_vector,Vector<IssmDouble> *vector){/*{{{*/ -
issm/trunk-jpl/src/c/classes/Nodes.cpp
r23642 r23926 4 4 */ 5 5 6 /*Headers : {{{*/6 /*Headers*/ 7 7 #ifdef HAVE_CONFIG_H 8 8 #include <config.h> … … 10 10 #error "Cannot compile with HAVE_CONFIG_H symbol! run configure first!" 11 11 #endif 12 13 12 #include "../shared/io/Comm/IssmComm.h" 14 13 #include "./Nodes.h" 15 14 #include "./Node.h" 16 17 15 using namespace std; 18 /*}}}*/19 16 20 17 /*Object constructors and destructor*/ … … 391 388 } 392 389 /*}}}*/ 390 391 void Nodes::CheckDofListAcrossPartitions(void){/*{{{*/ 392 393 /*recover my_rank:*/ 394 ISSM_MPI_Status status; 395 int my_rank = IssmComm::GetRank(); 396 int num_procs = IssmComm::GetSize(); 397 398 /*Display message*/ 399 if(VerboseModule()) _printf0_(" Checking degrees of freedom across partitions\n"); 400 401 /*Allocate vector to check degrees of freedom*/ 402 int gsize = this->NumberOfDofs(GsetEnum); 403 int glocalsize = this->NumberOfDofsLocal(GsetEnum); 404 Vector<IssmDouble>* dofs_check=new Vector<IssmDouble>(glocalsize,gsize); 405 406 /*First, go over all nodes, and masters can write their f dof and -1 for s-set*/ 407 for(int i=0;i<this->Size();i++){ 408 Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i)); 409 410 /*Skip if clone (will check later)*/ 411 if(node->IsClone()) continue; 412 413 /*Write degree of freedom if active*/ 414 int count = 0; 415 for(int j=0;j<node->gsize;j++){ 416 if(node->f_set[j]){ 417 if(node->s_set[j]) _error_("a degree of freedom is both in f and s set!"); 418 dofs_check->SetValue(node->gdoflist[j],reCast<IssmDouble>(node->fdoflist[count]),INS_VAL); 419 count++; 420 } 421 else{ 422 if(node->s_set[j]==0) _error_("a degree of freedom is neither in f nor in s set!"); 423 dofs_check->SetValue(node->gdoflist[j],-1.,INS_VAL); 424 } 425 } 426 } 427 dofs_check->Assemble(); 428 429 /*Get local vector with both masters and slaves:*/ 430 IssmDouble *local_dofs_check = NULL; 431 this->GetLocalVectorWithClonesGset(&local_dofs_check,dofs_check); 432 delete dofs_check; 433 434 /*Second, go over all nodes, and check that we still have what's expected...*/ 435 for(int i=0;i<this->Size();i++){ 436 Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i)); 437 438 /*Write degree of freedom if active*/ 439 int countg = 0; 440 int countf = 0; 441 int counts = 0; 442 for(int j=0;j<node->gsize;j++){ 443 int index = node->gdoflist_local[countg]; 444 if(node->f_set[j]){ 445 if(reCast<int>(local_dofs_check[index]) != node->fdoflist[countf]){ 446 _error_("Dof #"<<j<<" of node sid "<<node->Sid()<<" not consistent: "<<local_dofs_check[index]<<"!="<<node->fdoflist[countf]); 447 } 448 countf++; 449 } 450 else{ 451 if(local_dofs_check[index] != -1.){ 452 _error_("Dof #"<<j<<" of node sid "<<node->Sid()<<" not consistently in s set"); 453 } 454 counts++; 455 } 456 countg++; 457 } 458 } 459 460 /*cleanup and return*/ 461 xDelete<IssmDouble>(local_dofs_check); 462 }/*}}}*/ 463 void Nodes::GetLocalVectorWithClonesGset(IssmDouble** plocal_ug,Vector<IssmDouble> *ug){/*{{{*/ 464 465 /*recover my_rank:*/ 466 ISSM_MPI_Status status; 467 int my_rank = IssmComm::GetRank(); 468 int num_procs = IssmComm::GetSize(); 469 470 /*retrieve node info*/ 471 int glocalsize = this->NumberOfDofsLocalAll(GsetEnum); 472 int glocalsize_masters = this->NumberOfDofsLocal(GsetEnum); 473 int maxdofspernode = this->MaxNumDofs(GsetEnum); 474 475 /*Get local vector of ug*/ 476 int *indices_ug_masters = NULL; 477 IssmDouble *local_ug_masters = NULL; 478 ug->GetLocalVector(&local_ug_masters,&indices_ug_masters); 479 _assert_(glocalsize_masters==indices_ug_masters[glocalsize_masters-1] - indices_ug_masters[0]+1); 480 xDelete<int>(indices_ug_masters); 481 482 /*Now, extend vectors to account for clones (make vectors longer, for clones at the end)*/ 483 IssmDouble *local_ug = xNew<IssmDouble>(glocalsize); 484 xMemCpy<IssmDouble>(local_ug,local_ug_masters,glocalsize_masters); 485 xDelete<IssmDouble>(local_ug_masters); 486 487 /*Now send and receive ug for nodes on partition edge*/ 488 #ifdef _HAVE_AD_ 489 IssmDouble* buffer = xNew<IssmDouble>(this->Size()*maxdofspernode,"t"); //only one alloc, "t" is required by adolc 490 #else 491 IssmDouble* buffer = xNew<IssmDouble>(this->Size()*maxdofspernode); 492 #endif 493 for(int rank=0;rank<num_procs;rank++){ 494 if(this->common_send[rank]){ 495 int numids = this->common_send[rank]; 496 for(int i=0;i<numids;i++){ 497 int master_lid = this->common_send_ids[rank][i]; 498 Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(master_lid)); 499 _assert_(!node->IsClone()); 500 for(int j=0;j<node->gsize;j++) buffer[i*maxdofspernode+j]=local_ug[node->gdoflist_local[j]]; 501 } 502 ISSM_MPI_Send(buffer,numids*maxdofspernode,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm()); 503 } 504 } 505 for(int rank=0;rank<num_procs;rank++){ 506 if(this->common_recv[rank]){ 507 int numids = this->common_recv[rank]; 508 ISSM_MPI_Recv(buffer,numids*maxdofspernode,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status); 509 for(int i=0;i<numids;i++){ 510 int master_lid = this->common_recv_ids[rank][i]; 511 Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(master_lid)); 512 for(int j=0;j<node->gsize;j++) local_ug[node->gdoflist_local[j]] = buffer[i*maxdofspernode+j]; 513 } 514 } 515 } 516 xDelete<IssmDouble>(buffer); 517 518 /*Assign output pointer*/ 519 *plocal_ug = local_ug; 520 }/*}}}*/ -
issm/trunk-jpl/src/c/classes/Nodes.h
r23642 r23926 3 3 4 4 #include "../datastructures/datastructures.h" 5 #include "../toolkits/toolkits.h" 5 6 class Parameters; 6 7 class Elements; … … 47 48 int NumberOfNodesLocalAll(void); 48 49 bool RequiresDofReindexing(void); 50 void CheckDofListAcrossPartitions(void); 51 void GetLocalVectorWithClonesGset(IssmDouble** plocal_vector,Vector<IssmDouble> *vector); 49 52 }; 50 53
Note:
See TracChangeset
for help on using the changeset viewer.