Changeset 23600
- Timestamp:
- 01/04/19 21:15:03 (6 years ago)
- Location:
- issm/trunk-jpl/src/c
- Files:
-
- 4 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk-jpl/src/c/classes/Node.cpp
r23599 r23600 500 500 } 501 501 /*}}}*/ 502 void 502 void Node::SetApproximation(int in_approximation){/*{{{*/ 503 503 504 504 this->approximation = in_approximation; … … 749 749 } 750 750 /*}}}*/ 751 void Node::ShowTrueDofs(int* truedofs, int ncols,int setenum){/*{{{*/ 752 753 int j; 754 755 /*Are we a clone? : */ 756 if(indexing.clone) return; 751 void Node::ShowTrueDofs(int* truedofs,int setenum){/*{{{*/ 752 753 _assert_(!this->indexing.clone); 757 754 758 755 /*Ok, we are not a clone, just plug our dofs into truedofs: */ 759 756 switch(setenum){ 760 757 case GsetEnum: 761 for( j=0;j<this->indexing.gsize;j++) truedofs[ncols*sid+j]=indexing.gdoflist[j];758 for(int j=0;j<this->indexing.gsize;j++) truedofs[j]=indexing.gdoflist[j]; 762 759 break; 763 760 case FsetEnum: 764 for( j=0;j<this->indexing.fsize;j++) truedofs[ncols*sid+j]=indexing.fdoflist[j];761 for(int j=0;j<this->indexing.fsize;j++) truedofs[j]=indexing.fdoflist[j]; 765 762 break; 766 763 case SsetEnum: 767 for( j=0;j<this->indexing.ssize;j++) truedofs[ncols*sid+j]=indexing.sdoflist[j];764 for(int j=0;j<this->indexing.ssize;j++) truedofs[j]=indexing.sdoflist[j]; 768 765 break; 769 766 default: … … 773 770 } 774 771 /*}}}*/ 775 void Node::UpdateCloneDofs(int* alltruedofs,int ncols,int setenum){/*{{{*/ 776 777 int j; 778 779 /*If we are not a clone, don't update, we already have dofs!: */ 780 if(!indexing.clone)return; 772 void Node::UpdateCloneDofs(int* alltruedofs,int setenum){/*{{{*/ 773 774 _assert_(this->indexing.clone); 781 775 782 776 /*Ok, we are a clone node, but we did not create the dofs for this node. … … 784 778 switch(setenum){ 785 779 case GsetEnum: 786 for( j=0;j<this->indexing.gsize;j++) indexing.gdoflist[j]=alltruedofs[ncols*sid+j];780 for(int j=0;j<this->indexing.gsize;j++) indexing.gdoflist[j]=alltruedofs[j]; 787 781 break; 788 782 case FsetEnum: 789 for( j=0;j<this->indexing.fsize;j++) indexing.fdoflist[j]=alltruedofs[ncols*sid+j];783 for(int j=0;j<this->indexing.fsize;j++) indexing.fdoflist[j]=alltruedofs[j]; 790 784 break; 791 785 case SsetEnum: 792 for( j=0;j<this->indexing.ssize;j++) indexing.sdoflist[j]=alltruedofs[ncols*sid+j];786 for(int j=0;j<this->indexing.ssize;j++) indexing.sdoflist[j]=alltruedofs[j]; 793 787 break; 794 788 default: -
issm/trunk-jpl/src/c/classes/Node.h
r23599 r23600 78 78 bool RequiresDofReindexing(void); 79 79 void SetCurrentConfiguration(DataSet* nodes,Vertices* vertices); 80 void ShowTrueDofs(int* truerows,int ncols,intsetenum);80 void ShowTrueDofs(int* truerows,int setenum); 81 81 int Sid(void); 82 void UpdateCloneDofs(int* alltruerows,int ncols,intsetenum);82 void UpdateCloneDofs(int* alltruerows,int setenum); 83 83 void VecMerge(Vector<IssmDouble>* ug, IssmDouble* vector_serial,int setenum); 84 84 void VecReduce(Vector<IssmDouble>* vector, IssmDouble* ug_serial,int setnum); -
issm/trunk-jpl/src/c/classes/Nodes.cpp
r23599 r23600 133 133 void Nodes::DistributeDofs(int setenum){/*{{{*/ 134 134 135 /*some check: */136 _assert_(setenum==GsetEnum || setenum==FsetEnum || setenum==SsetEnum);137 138 int i;139 int dofcount=0;140 int maxdofspernode=0;141 int* alldofcount=NULL;142 int* truedofs=NULL;143 int* alltruedofs=NULL;144 int numnodes=0;145 146 135 /*recover my_rank:*/ 136 ISSM_MPI_Status status; 147 137 int my_rank = IssmComm::GetRank(); 148 138 int num_procs = IssmComm::GetSize(); 149 139 150 140 /*Go through objects, and distribute dofs locally, from 0 to numberofdofsperobject*/ 151 for(i=0;i<this->Size();i++){ 141 int dofcount=0; 142 for(int i=0;i<this->Size();i++){ 152 143 Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i)); 153 144 node->DistributeDofs(&dofcount,setenum); … … 158 149 * cpus by the total last dofs of the previus cpu, starting from 0. 159 150 * First: get number of dofs for each cpu*/ 160 alldofcount=xNew<int>(num_procs);151 int* alldofcount=xNew<int>(num_procs); 161 152 ISSM_MPI_Gather(&dofcount,1,ISSM_MPI_INT,alldofcount,1,ISSM_MPI_INT,0,IssmComm::GetComm()); 162 153 ISSM_MPI_Bcast(alldofcount,num_procs,ISSM_MPI_INT,0,IssmComm::GetComm()); 163 154 164 155 /* Every cpu should start its own dof count at the end of the dofcount from cpu-1*/ 165 dofcount=0;166 for(i =0;i<my_rank;i++){167 dofcount+=alldofcount[i];168 } 169 for(i =0;i<this->Size();i++){156 int offset=0; 157 for(int i=0;i<my_rank;i++) offset+=alldofcount[i]; 158 xDelete<int>(alldofcount); 159 160 for(int i=0;i<this->Size();i++){ 170 161 /*Check that this node corresponds to our analysis currently being carried out: */ 171 162 Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i)); 172 node->OffsetDofs( dofcount,setenum);163 node->OffsetDofs(offset,setenum); 173 164 } 174 165 … … 176 167 * object that is not a clone, tell them to show their dofs, so that later on, they can get picked 177 168 * up by their clones: */ 178 maxdofspernode=this->MaxNumDofs(setenum); 179 numnodes=this->NumberOfNodes(); 180 if(numnodes*maxdofspernode){ 181 truedofs= xNewZeroInit<int>(numnodes*maxdofspernode); //initialize to 0, so that we can pick up the max 182 alltruedofs=xNewZeroInit<int>(numnodes*maxdofspernode); 183 } 184 185 for(i=0;i<this->Size();i++){ 186 Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i)); 187 node->ShowTrueDofs(truedofs,maxdofspernode,setenum);//give maxdofspernode, column size, so that nodes can index into truedofs 188 } 189 190 ISSM_MPI_Allreduce((void*)truedofs,(void*)alltruedofs,numnodes*maxdofspernode,ISSM_MPI_INT,ISSM_MPI_MAX,IssmComm::GetComm()); 191 192 /* Now every cpu knows the true dofs of everyone else that is not a clone*/ 193 for(i=0;i<this->Size();i++){ 194 Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i)); 195 node->UpdateCloneDofs(alltruedofs,maxdofspernode,setenum); 196 } 169 int maxdofspernode = this->MaxNumDofs(setenum); 170 int* truedofs = xNew<int>(this->Size()*maxdofspernode); //only one alloc 171 for(int rank=0;rank<num_procs;rank++){ 172 if(this->common_send[rank]){ 173 int numids = this->common_send[rank]; 174 for(int i=0;i<numids;i++){ 175 Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(this->common_send_ids[rank][i])); 176 node->ShowTrueDofs(&truedofs[i*maxdofspernode+0],setenum); 177 } 178 ISSM_MPI_Send(truedofs,numids*maxdofspernode,ISSM_MPI_INT,rank,0,IssmComm::GetComm()); 179 } 180 } 181 for(int rank=0;rank<num_procs;rank++){ 182 if(this->common_recv[rank]){ 183 int numids = this->common_recv[rank]; 184 ISSM_MPI_Recv(truedofs,numids*maxdofspernode,ISSM_MPI_INT,rank,0,IssmComm::GetComm(),&status); 185 for(int i=0;i<numids;i++){ 186 Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(this->common_recv_ids[rank][i])); 187 node->UpdateCloneDofs(&truedofs[i*maxdofspernode+0],setenum); 188 } 189 } 190 } 191 xDelete<int>(truedofs); 197 192 198 193 /*Update indexingupdateflag*/ 199 for(i =0;i<this->Size();i++){194 for(int i=0;i<this->Size();i++){ 200 195 Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i)); 201 196 node->ReindexingDone(); 202 197 } 203 204 /* Free ressources: */205 xDelete<int>(alldofcount);206 xDelete<int>(truedofs);207 xDelete<int>(alltruedofs);208 198 } 209 199 /*}}}*/ -
issm/trunk-jpl/src/c/modules/NodesDofx/NodesDofx.cpp
r23587 r23600 24 24 nodes->DistributeDofs(FsetEnum); 25 25 nodes->DistributeDofs(SsetEnum); 26 27 26 }
Note:
See TracChangeset
for help on using the changeset viewer.