Changeset 23638
- Timestamp:
- 01/16/19 16:54:18 (6 years ago)
- Location:
- issm/trunk-jpl/src/c
- Files:
-
- 10 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk-jpl/src/c/classes/Elements/Penta.cpp
r23629 r23638 1662 1662 1663 1663 switch(type){ 1664 case VertexLIdEnum: 1665 for (int i=0;i<NUMVERTICES;i++){ 1666 values[i]=vector[this->vertices[i]->Lid()]; 1667 } 1668 /*update input*/ 1669 this->inputs->AddInput(new PentaInput(name,values,P1Enum)); 1670 return; 1671 1664 1672 case VertexPIdEnum: 1665 1673 for (int i=0;i<NUMVERTICES;i++){ -
issm/trunk-jpl/src/c/classes/FemModel.cpp
r23636 r23638 1364 1364 *plocal_ug = local_ug; 1365 1365 }/*}}}*/ 1366 void FemModel::GetLocalVectorWithClonesVertices(IssmDouble** plocal_vector,Vector<IssmDouble> *vector){/*{{{*/ 1367 1368 /*recover my_rank:*/ 1369 ISSM_MPI_Status status; 1370 int my_rank = IssmComm::GetRank(); 1371 int num_procs = IssmComm::GetSize(); 1372 1373 /*retrieve vertex info*/ 1374 int localsize = this->vertices->NumberOfVerticesLocalAll(); 1375 int localsize_masters = this->vertices->NumberOfVerticesLocal(); 1376 1377 /*Get local vector of vector*/ 1378 int *indices_vector_masters = NULL; 1379 IssmDouble *local_vector_masters = NULL; 1380 vector->GetLocalVector(&local_vector_masters,&indices_vector_masters); 1381 _assert_(localsize_masters==indices_vector_masters[localsize_masters-1] - indices_vector_masters[0]+1); 1382 xDelete<int>(indices_vector_masters); 1383 1384 /*Now, extend vectors to account for clones (make vectors longer, for clones at the end)*/ 1385 IssmDouble *local_vector = xNew<IssmDouble>(localsize); 1386 xMemCpy<IssmDouble>(local_vector,local_vector_masters,localsize_masters); 1387 xDelete<IssmDouble>(local_vector_masters); 1388 1389 /*Now send and receive vector for vertices on partition edge*/ 1390 #ifdef _HAVE_AD_ 1391 IssmDouble* buffer = xNew<IssmDouble>(this->vertices->Size(),"t"); //only one alloc, "t" is required by adolc 1392 #else 1393 IssmDouble* buffer = xNew<IssmDouble>(this->vertices->Size()); 1394 #endif 1395 for(int rank=0;rank<num_procs;rank++){ 1396 if(this->vertices->common_send[rank]){ 1397 int numids = this->vertices->common_send[rank]; 1398 for(int i=0;i<numids;i++){ 1399 int master_lid = this->vertices->common_send_ids[rank][i]; 1400 Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid)); 1401 _assert_(!vertex->clone); 1402 buffer[i] = local_vector[vertex->lid]; 1403 } 1404 ISSM_MPI_Send(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm()); 1405 } 1406 } 1407 for(int rank=0;rank<num_procs;rank++){ 1408 if(this->vertices->common_recv[rank]){ 1409 int numids = this->vertices->common_recv[rank]; 1410 ISSM_MPI_Recv(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status); 1411 for(int i=0;i<numids;i++){ 1412 int master_lid = this->vertices->common_recv_ids[rank][i]; 1413 Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid)); 1414 _assert_(vertex->clone); 1415 local_vector[vertex->lid] = buffer[i]; 1416 } 1417 } 1418 } 1419 xDelete<IssmDouble>(buffer); 1420 1421 /*Assign output pointer*/ 1422 *plocal_vector = local_vector; 1423 }/*}}}*/ 1366 1424 void FemModel::GroundedAreax(IssmDouble* pV, bool scaled){/*{{{*/ 1367 1425 … … 2616 2674 2617 2675 /*Allocate vector*/ 2618 Vector<IssmDouble> *vx=new Vector<IssmDouble>(vertices->NumberOfVertices()); 2619 Vector<IssmDouble> *vy=new Vector<IssmDouble>(vertices->NumberOfVertices()); 2620 Vector<IssmDouble> *vz=new Vector<IssmDouble>(vertices->NumberOfVertices()); 2676 int numvert = vertices->NumberOfVertices(); 2677 int numvert_local = vertices->NumberOfVerticesLocal(); 2678 Vector<IssmDouble> *vx=new Vector<IssmDouble>(numvert_local,numvert); 2679 Vector<IssmDouble> *vy=new Vector<IssmDouble>(numvert_local,numvert); 2680 Vector<IssmDouble> *vz=new Vector<IssmDouble>(numvert_local,numvert); 2621 2681 2622 2682 /*Update verices new geometry: */ -
issm/trunk-jpl/src/c/classes/FemModel.h
r23636 r23638 97 97 void GetInputLocalMinMaxOnNodesx(IssmDouble** pmin,IssmDouble** pmax,IssmDouble* ug); 98 98 void GetLocalVectorWithClonesGset(IssmDouble** plocal_ug,Vector<IssmDouble> *ug); 99 void GetLocalVectorWithClonesVertices(IssmDouble** plocal_vector,Vector<IssmDouble> *vector); 99 100 void GroundedAreax(IssmDouble* pV, bool scaled); 100 101 void IceMassx(IssmDouble* pV, bool scaled); -
issm/trunk-jpl/src/c/classes/Nodes.cpp
r23629 r23638 171 171 172 172 for(int i=0;i<this->Size();i++){ 173 /*Check that this node corresponds to our analysis currently being carried out: */174 173 Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i)); 175 174 node->DistributeGlobalDofsMasters(offset,setenum); -
issm/trunk-jpl/src/c/classes/Vertex.cpp
r23599 r23638 20 20 } 21 21 /*}}}*/ 22 Vertex::Vertex(int vertex_id, int vertex_sid,int vertex_ lid,int vertex_pid,bool vertex_clone, IoModel* iomodel,bool isamr){/*{{{*/22 Vertex::Vertex(int vertex_id, int vertex_sid,int vertex_pid,bool vertex_clone, IoModel* iomodel,bool isamr){/*{{{*/ 23 23 24 24 /*Checks in debugging mode*/ … … 29 29 this->sid = vertex_sid; 30 30 this->pid = vertex_pid; 31 this->lid = vertex_lid;31 this->lid = -1; /*Assigned later*/ 32 32 this->clone = vertex_clone; 33 33 -
issm/trunk-jpl/src/c/classes/Vertex.h
r23599 r23638 37 37 /*Vertex constructors, destructors {{{*/ 38 38 Vertex(); 39 Vertex(int id, int sid,int lid,intpid,bool clone, IoModel* iomodel,bool isamr);39 Vertex(int id, int sid,int pid,bool clone, IoModel* iomodel,bool isamr); 40 40 ~Vertex(); 41 41 /*}}}*/ -
issm/trunk-jpl/src/c/classes/Vertices.cpp
r23599 r23638 25 25 /*Object constructors and destructor*/ 26 26 Vertices::Vertices(){/*{{{*/ 27 this->enum_type = VerticesEnum; 28 this->common_recv = NULL; 29 this->common_recv_ids = NULL; 30 this->common_send = NULL; 31 this->common_send_ids = NULL; 27 this->enum_type = VerticesEnum; 28 this->common_recv = NULL; 29 this->common_recv_ids = NULL; 30 this->common_send = NULL; 31 this->common_send_ids = NULL; 32 this->numberofvertices = -1; 33 this->numberofvertices_local = -1; 34 this->numberofmasters_local = -1; 32 35 return; 33 36 } … … 127 130 } 128 131 /*}}}*/ 129 int Vertices::NumberOfVertices(void){/*{{{*/130 131 int i,sid;132 int max_sid=0;133 int vertex_max_sid;134 135 if(this->Size()==0) return 0;136 137 for(i=0;i<this->Size();i++){138 Vertex* vertex=xDynamicCast<Vertex*>(this->GetObjectByOffset(i));139 sid=vertex->Sid();140 if (sid>max_sid)max_sid=sid;141 }142 143 ISSM_MPI_Reduce (&max_sid,&vertex_max_sid,1,ISSM_MPI_INT,ISSM_MPI_MAX,0,IssmComm::GetComm() );144 ISSM_MPI_Bcast(&vertex_max_sid,1,ISSM_MPI_INT,0,IssmComm::GetComm());145 max_sid=vertex_max_sid;146 147 /*sid starts at 0*/148 max_sid++;149 150 /*return:*/151 return max_sid;152 }153 /*}}}*/154 132 void Vertices::LatLonList(IssmDouble** plat,IssmDouble** plon){/*{{{*/ 155 133 … … 189 167 } 190 168 /*}}}*/ 169 170 void Vertices::Finalize(){/*{{{*/ 171 172 /*recover my_rank:*/ 173 ISSM_MPI_Status status; 174 int my_rank = IssmComm::GetRank(); 175 int num_procs = IssmComm::GetSize(); 176 177 /*1. set counters*/ 178 this->numberofvertices_local=this->Size(); 179 this->numberofmasters_local=0; 180 for(int i=0;i<this->Size();i++){ 181 Vertex* vertex=xDynamicCast<Vertex*>(this->GetObjectByOffset(i)); 182 if(!vertex->clone) this->numberofmasters_local++; 183 } 184 ISSM_MPI_Allreduce((void*)&this->numberofmasters_local,(void*)&this->numberofvertices,1,ISSM_MPI_INT,ISSM_MPI_SUM,IssmComm::GetComm()); 185 186 /*2. Distribute lids (First: masters, then clones)*/ 187 int lid = 0; 188 for(int i=0;i<this->Size();i++){ 189 Vertex* vertex=xDynamicCast<Vertex*>(this->GetObjectByOffset(i)); 190 if(!vertex->clone) vertex->lid=lid++; 191 } 192 for(int i=0;i<this->Size();i++){ 193 Vertex* vertex=xDynamicCast<Vertex*>(this->GetObjectByOffset(i)); 194 if(vertex->clone) vertex->lid=lid++; 195 } 196 197 /* Now every object has distributed dofs, but locally, and with a dof count starting from 198 * 0. This means the dofs between all the cpus are not unique. We now offset the dofs of each 199 * cpus by the total last (master) dofs of the previus cpu, starting from 0. 200 * First: get number of dofs for each cpu*/ 201 int* all_num_masters=xNew<int>(num_procs); 202 ISSM_MPI_Gather(&this->numberofmasters_local,1,ISSM_MPI_INT,all_num_masters,1,ISSM_MPI_INT,0,IssmComm::GetComm()); 203 ISSM_MPI_Bcast(all_num_masters,num_procs,ISSM_MPI_INT,0,IssmComm::GetComm()); 204 205 /* Every cpu should start its own dof count at the end of the dofcount from cpu-1*/ 206 int offset=0; 207 for(int i=0;i<my_rank;i++) offset+=all_num_masters[i]; 208 xDelete<int>(all_num_masters); 209 for(int i=0;i<this->Size();i++){ 210 Vertex* vertex=xDynamicCast<Vertex*>(this->GetObjectByOffset(i)); 211 vertex->pid = vertex->lid+offset; 212 } 213 214 /* Finally, remember that cpus may have skipped some objects, because they were clones. For every 215 * object that is not a clone, tell them to show their dofs, so that later on, they can get picked 216 * up by their clones: */ 217 int* truepids = xNew<int>(this->Size()); //only one alloc 218 for(int rank=0;rank<num_procs;rank++){ 219 if(this->common_send[rank]){ 220 int numids = this->common_send[rank]; 221 for(int i=0;i<numids;i++){ 222 Vertex* vertex=xDynamicCast<Vertex*>(this->GetObjectByOffset(this->common_send_ids[rank][i])); 223 truepids[i] = vertex->pid; 224 } 225 ISSM_MPI_Send(truepids,numids,ISSM_MPI_INT,rank,0,IssmComm::GetComm()); 226 } 227 } 228 for(int rank=0;rank<num_procs;rank++){ 229 if(this->common_recv[rank]){ 230 int numids = this->common_recv[rank]; 231 ISSM_MPI_Recv(truepids,numids,ISSM_MPI_INT,rank,0,IssmComm::GetComm(),&status); 232 for(int i=0;i<numids;i++){ 233 Vertex* vertex=xDynamicCast<Vertex*>(this->GetObjectByOffset(this->common_recv_ids[rank][i])); 234 vertex->pid = truepids[i]; 235 } 236 } 237 } 238 xDelete<int>(truepids); 239 }/*}}}*/ 240 int Vertices::NumberOfVertices(){/*{{{*/ 241 return this->numberofvertices; 242 }/*}}}*/ 243 int Vertices::NumberOfVerticesLocal(void){/*{{{*/ 244 return this->numberofmasters_local; 245 }/*}}}*/ 246 int Vertices::NumberOfVerticesLocalAll(void){/*{{{*/ 247 return this->numberofvertices_local; 248 }/*}}}*/ -
issm/trunk-jpl/src/c/classes/Vertices.h
r23599 r23638 14 14 class Vertices: public DataSet{ 15 15 16 private: 17 int numberofvertices; 18 int numberofvertices_local; 19 int numberofmasters_local; 16 20 public: 17 21 int* common_recv; … … 26 30 /*Objects virtual functions*/ 27 31 Vertices* Copy(); 28 void Marshall(char** pmarshalled_data,int* pmarshalled_data_size, int marshall_direction);32 void Marshall(char** pmarshalled_data,int* pmarshalled_data_size, int marshall_direction); 29 33 30 34 /*numerics:*/ 35 void Finalize(void); 31 36 int NumberOfVertices(void); 37 int NumberOfVerticesLocal(void); 38 int NumberOfVerticesLocalAll(void); 32 39 void LatLonList(IssmDouble** lat,IssmDouble** lon); 33 40 }; -
issm/trunk-jpl/src/c/modules/InputUpdateFromVectorx/InputUpdateFromVectorx.cpp
r23524 r23638 9 9 void InputUpdateFromVectorx(FemModel* femmodel,Vector<IssmDouble>* vector, int name, int type){ 10 10 11 IssmDouble* serial_vector=vector->ToMPISerial(); 12 InputUpdateFromVectorx(femmodel,serial_vector,name,type); 13 xDelete<IssmDouble>(serial_vector); 11 if(type==VertexPIdEnum){ 12 IssmDouble* serial_vector=NULL; 13 femmodel->GetLocalVectorWithClonesVertices(&serial_vector,vector); 14 InputUpdateFromVectorx(femmodel,serial_vector,name,VertexLIdEnum); 15 xDelete<IssmDouble>(serial_vector); 16 } 17 else{ 18 IssmDouble* serial_vector=vector->ToMPISerial(); 19 InputUpdateFromVectorx(femmodel,serial_vector,name,type); 20 xDelete<IssmDouble>(serial_vector); 21 } 14 22 } 15 23 -
issm/trunk-jpl/src/c/modules/ModelProcessorx/CreateElementsVerticesAndMaterials.cpp
r23635 r23638 411 411 if(vertices_lids[i]!=-1){ 412 412 bool isclone = (vertices_ranks[MAXCONNECTIVITY*i+0]!=my_rank); 413 vertices->AddObject(new Vertex(i+1,i,vertices_ lids[i],vertices_pids[i],isclone,iomodel,isamr));413 vertices->AddObject(new Vertex(i+1,i,vertices_pids[i],isclone,iomodel,isamr)); 414 414 } 415 415 } … … 424 424 if(vertices_lids[i]!=-1){ 425 425 bool isclone = (vertices_ranks[MAXCONNECTIVITY*i+0]!=my_rank); 426 vertices->AddObject(new Vertex(i+1,i,vertices_ lids[i],vertices_pids[i],isclone,iomodel,isamr));426 vertices->AddObject(new Vertex(i+1,i,vertices_pids[i],isclone,iomodel,isamr)); 427 427 } 428 428 } … … 453 453 vertices->common_send_ids=common_send_ids; 454 454 vertices->common_recv_ids=common_recv_ids; 455 }/*}}}*/ 455 456 /*Finalize Initialization*/ 457 vertices->Finalize(); 458 }/*}}}*/
Note:
See TracChangeset
for help on using the changeset viewer.