Changeset 15861
- Timestamp:
- 08/21/13 15:44:01 (12 years ago)
- Location:
- issm/trunk-jpl/src/c
- Files:
-
- 5 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk-jpl/src/c/classes/Node.cpp
r15771 r15861 39 39 40 40 /*indexing:*/ 41 this->indexingupdate = true; 41 42 DistributeNumDofs(&this->indexing,analysis_type,in_approximation); //number of dofs per node 42 43 gsize=this->indexing.gsize; … … 127 128 _printf_(" analysis_type: " << EnumToStringx(analysis_type) << "\n"); 128 129 _printf_(" approximation: " << EnumToStringx(approximation) << "\n"); 130 _printf_(" indexingupdate: " << indexingupdate << "\n"); 129 131 indexing.Echo(); 130 132 _printf_(" inputs: " << inputs << "\n"); … … 139 141 _printf_(" sid: " << sid << "\n"); 140 142 _printf_(" analysis_type: " << EnumToStringx(analysis_type) << "\n"); 143 _printf_(" approximation: " << EnumToStringx(approximation) << "\n"); 144 _printf_(" indexingupdate: " << indexingupdate << "\n"); 141 145 indexing.DeepEcho(); 142 146 _printf_(" inputs\n"); … … 159 163 int Node::GetDof(int dofindex,int setenum){ 160 164 165 _assert_(!this->indexingupdate); 161 166 if(setenum==GsetEnum){ 162 167 _assert_(dofindex>=0 && dofindex<indexing.gsize); … … 179 184 int count=0; 180 185 int count2=0; 186 187 _assert_(!this->indexingupdate); 181 188 182 189 if(approximation_enum==NoneApproximationEnum){ … … 242 249 int count2=0; 243 250 251 _assert_(!this->indexingupdate); 252 244 253 if(approximation_enum==NoneApproximationEnum){ 245 254 if(setenum==GsetEnum)for(i=0;i<this->indexing.gsize;i++) outdoflist[i]=i; … … 365 374 } 366 375 /*}}}*/ 376 /*FUNCTION Node::RequiresDofReindexing{{{*/ 377 bool Node::RequiresDofReindexing(void){ 378 379 return this->indexingupdate; 380 381 } 382 /*}}}*/ 383 /*FUNCTION Node::ReindexingDone{{{*/ 384 void Node::ReindexingDone(void){ 385 386 this->indexingupdate = false; 387 388 } 389 /*}}}*/ 367 390 /*FUNCTION Node::RelaxConstraint{{{*/ 368 391 void Node::RelaxConstraint(int dof){ … … 376 399 void Node::CreateVecSets(Vector<IssmDouble>* pv_g,Vector<IssmDouble>* pv_f,Vector<IssmDouble>* pv_s){ 377 400 378 IssmDouble gvalue=1. 0; //all nodes are in the g set;401 IssmDouble gvalue=1.; //all nodes are in the g set; 379 402 IssmDouble value; 380 403 381 int i; 382 383 for(i=0;i<this->indexing.gsize;i++){ 404 for(int i=0;i<this->indexing.gsize;i++){ 384 405 385 406 /*g set: */ … … 457 478 void Node::Deactivate(void){ 458 479 459 indexing.Deactivate(); 480 if(IsActive()){ 481 this->indexingupdate = true; 482 indexing.Deactivate(); 483 } 460 484 461 485 } … … 464 488 void Node::Activate(void){ 465 489 466 indexing.Activate(); 490 if(!IsActive()){ 491 this->indexingupdate = true; 492 indexing.Activate(); 493 } 467 494 468 495 } … … 491 518 492 519 if(approximation_enum==NoneApproximationEnum){ 493 if (setenum==GsetEnum) numdofs=this->indexing.gsize;520 if (setenum==GsetEnum) numdofs=this->indexing.gsize; 494 521 else if (setenum==FsetEnum) numdofs=this->indexing.fsize; 495 522 else if (setenum==SsetEnum) numdofs=this->indexing.ssize; … … 621 648 void Node::UpdateSpcs(IssmDouble* ys){ 622 649 623 int count=0; 624 int i; 625 626 count=0; 627 for(i=0;i<this->indexing.gsize;i++){ 650 int count = 0; 651 for(int i=0;i<this->indexing.gsize;i++){ 628 652 if(this->indexing.s_set[i]){ 629 653 this->indexing.svalues[i]=ys[this->indexing.sdoflist[count]]; … … 634 658 /*}}}*/ 635 659 /*FUNCTION Node::VecMerge {{{*/ 636 void 637 638 IssmDouble * values=NULL;639 int * indices=NULL;640 int count=0;641 int i;660 void Node::VecMerge(Vector<IssmDouble>* ug, IssmDouble* vector_serial,int setenum){ 661 662 IssmDouble *values = NULL; 663 int *indices = NULL; 664 int count = 0; 665 int i; 642 666 643 667 if(setenum==FsetEnum){ -
issm/trunk-jpl/src/c/classes/Node.h
r15771 r15861 33 33 int sid; //"serial" id (rank of this node if the dataset was serial on 1 cpu) 34 34 35 bool indexingupdate; 35 36 DofIndexing indexing; 36 37 Inputs *inputs; //properties of this node … … 88 89 void Deactivate(void); 89 90 void UpdateSpcs(IssmDouble* ys); 91 void ReindexingDone(void); 92 bool RequiresDofReindexing(void); 90 93 int IsFloating(); 91 94 int IsGrounded(); -
issm/trunk-jpl/src/c/classes/Nodes.cpp
r15838 r15861 51 51 52 52 /*Go through objects, and distribute dofs locally, from 0 to numberofdofsperobject*/ 53 for 54 Node* node=dynamic_cast<Node*>(this->GetObjectByOffset(i)); 55 56 /*Check that this node corresponds to our analysis currently being carried out: */ 57 if 53 for(i=0;i<this->Size();i++){ 54 Node* node=dynamic_cast<Node*>(this->GetObjectByOffset(i)); 55 56 /*Check that this node corresponds to our analysis currently being carried out: */ 57 if(node->InAnalysis(analysis_type)){ 58 58 node->DistributeDofs(&dofcount,setenum); 59 59 } … … 73 73 dofcount+=alldofcount[i]; 74 74 } 75 for 75 for(i=0;i<this->Size();i++){ 76 76 /*Check that this node corresponds to our analysis currently being carried out: */ 77 77 Node* node=dynamic_cast<Node*>(this->GetObjectByOffset(i)); … … 91 91 } 92 92 93 for 93 for(i=0;i<this->Size();i++){ 94 94 Node* node=dynamic_cast<Node*>(this->GetObjectByOffset(i)); 95 95 if (node->InAnalysis(analysis_type)){ … … 101 101 102 102 /* Now every cpu knows the true dofs of everyone else that is not a clone*/ 103 for 103 for(i=0;i<this->Size();i++){ 104 104 Node* node=dynamic_cast<Node*>(this->GetObjectByOffset(i)); 105 105 if (node->InAnalysis(analysis_type)){ 106 106 node->UpdateCloneDofs(alltruedofs,maxdofspernode,setenum); 107 } 108 } 109 110 /*Update indexingupdateflag*/ 111 for(i=0;i<this->Size();i++){ 112 Node* node=dynamic_cast<Node*>(this->GetObjectByOffset(i)); 113 if (node->InAnalysis(analysis_type)){ 114 node->ReindexingDone(); 107 115 } 108 116 } … … 330 338 } 331 339 /*}}}*/ 340 /*FUNCTION Nodes::RequiresDofReindexing{{{*/ 341 bool Nodes::RequiresDofReindexing(int analysis_type){ 342 343 int flag = 0; 344 int allflag; 345 346 /*Now go through all nodes, and get how many dofs they own, unless they are clone nodes: */ 347 for(int i=0;i<this->Size();i++){ 348 349 Node* node=dynamic_cast<Node*>(this->GetObjectByOffset(i)); 350 351 /*Check that this node corresponds to our analysis currently being carried out: */ 352 if(node->InAnalysis(analysis_type)){ 353 if(node->RequiresDofReindexing()){ 354 flag = 1; 355 break; 356 } 357 } 358 } 359 360 /*Grab max of all cpus: */ 361 ISSM_MPI_Allreduce((void*)&flag,(void*)&allflag,1,ISSM_MPI_INT,ISSM_MPI_MAX,IssmComm::GetComm()); 362 363 if(allflag){ 364 return true; 365 } 366 else{ 367 return false; 368 } 369 } 370 /*}}}*/ -
issm/trunk-jpl/src/c/classes/Nodes.h
r15634 r15861 27 27 void DistributeDofs(int analysis_type,int SETENUM); 28 28 void FlagClones(int analysis_type); 29 bool RequiresDofReindexing(int analysis_type); 29 30 int MaxNumDofs(int analysis_type,int setenum); 30 31 int MaximumId(void); -
issm/trunk-jpl/src/c/modules/NodesDofx/NodesDofx.cpp
r14999 r15861 10 10 void NodesDofx(Nodes* nodes, Parameters* parameters,int configuration_type){ 11 11 12 int noerr=1;13 i nt found=0;12 /*Do we have any nodes for this analysis type? :*/ 13 if(!nodes->NumberOfNodes(configuration_type)) return; 14 14 15 /*Do we have any nodes for this analysis type? :*/16 if( nodes->NumberOfNodes(configuration_type)){15 /*Do we really need to update dof indexings*/ 16 if(!nodes->RequiresDofReindexing(configuration_type)) return; 17 17 18 /*Ensure that only for each cpu, the partition border nodes only will be taken into account once 19 * across the cluster. To do so, we flag all the clone nodes: */ 20 nodes->FlagClones(configuration_type); 18 if(VerboseModule()) _printf0_(" Renumbering degrees of freedom\n"); 21 19 22 /*Go through all nodes, and build degree of freedom lists. Each node gets a fixed number of dofs. When 23 *a node has already been distributed dofs on one cpu, all other cpus with the same node cannot distribute it 24 *anymore. Use clone field to be sure of that: */ 25 nodes->DistributeDofs(configuration_type,GsetEnum); 26 nodes->DistributeDofs(configuration_type,FsetEnum); 27 nodes->DistributeDofs(configuration_type,SsetEnum); 28 } 20 /*Ensure that only for each cpu, the partition border nodes only will be taken into account once 21 * across the cluster. To do so, we flag all the clone nodes: */ 22 nodes->FlagClones(configuration_type); 23 24 /*Go through all nodes, and build degree of freedom lists. Each node gets a fixed number of dofs. When 25 *a node has already been distributed dofs on one cpu, all other cpus with the same node cannot distribute it 26 *anymore. Use clone field to be sure of that: */ 27 nodes->DistributeDofs(configuration_type,GsetEnum); 28 nodes->DistributeDofs(configuration_type,FsetEnum); 29 nodes->DistributeDofs(configuration_type,SsetEnum); 29 30 30 31 }
Note:
See TracChangeset
for help on using the changeset viewer.