Ice Sheet System Model  4.18
Code documentation
Public Member Functions | Data Fields | Private Attributes
Vertices Class Reference

Declaration of Vertices class. More...

#include <Vertices.h>

Inheritance diagram for Vertices:
DataSet

Public Member Functions

 Vertices ()
 
 ~Vertices ()
 
VerticesCopy ()
 
void Marshall (char **pmarshalled_data, int *pmarshalled_data_size, int marshall_direction)
 
void Finalize (IoModel *iomodel)
 
int NumberOfVertices (void)
 
int NumberOfVerticesLocal (void)
 
int NumberOfVerticesLocalAll (void)
 
void LatLonList (IssmDouble **lat, IssmDouble **lon)
 
- Public Member Functions inherited from DataSet
 DataSet ()
 
 DataSet (int enum_type)
 
 ~DataSet ()
 
void Marshall (char **pmarshalled_data, int *pmarshalled_data_size, int marshall_direction)
 
int GetEnum ()
 
int GetEnum (int offset)
 
void Echo ()
 
void DeepEcho ()
 
int AddObject (Object *object)
 
int DeleteObject (int id)
 
int Size ()
 
void clear ()
 
ObjectGetObjectByOffset (int offset)
 
ObjectGetObjectById (int *poffset, int eid)
 
void Presort ()
 
void Sort ()
 
DataSetCopy (void)
 
int DeleteObject (Object *object)
 

Data Fields

int * common_recv
 
int ** common_recv_ids
 
int * common_send
 
int ** common_send_ids
 
- Data Fields inherited from DataSet
std::vector< Object * > objects
 
int enum_type
 
int sorted
 
int presorted
 
int numsorted
 
int * sorted_ids
 
int * id_offsets
 

Private Attributes

int numberofvertices
 
int numberofvertices_local
 
int numberofmasters_local
 

Detailed Description

Declaration of Vertices class.

Declaration of Vertices class. Vertices are vector lists (Containers) of Vertex objects. A vertex is a set of (x,y,z) coordinates defining the location of points in the mesh (not to be confused with a node, which is a degree of freedom (DOF) for a particular analysis).

Definition at line 15 of file Vertices.h.

Constructor & Destructor Documentation

◆ Vertices()

Vertices::Vertices ( )

Definition at line 26 of file Vertices.cpp.

26  {/*{{{*/
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;
35  return;
36 }

◆ ~Vertices()

Vertices::~Vertices ( )

Definition at line 38 of file Vertices.cpp.

38  {/*{{{*/
39 
40  int num_proc=IssmComm::GetSize();
41 
42  if(this->common_recv) xDelete<int>(common_recv);
43  if(this->common_send) xDelete<int>(common_send);
44  if(this->common_recv_ids){
45  for(int i=0;i<num_proc;i++) if(common_recv_ids[i]) xDelete<int>(common_recv_ids[i]);
46  xDelete<int*>(common_recv_ids);
47  }
48  if(this->common_send_ids){
49  for(int i=0;i<num_proc;i++) if(common_send_ids[i]) xDelete<int>(common_send_ids[i]);
50  xDelete<int*>(common_send_ids);
51  }
52  return;
53 }

Member Function Documentation

◆ Copy()

Vertices * Vertices::Copy ( void  )

Definition at line 57 of file Vertices.cpp.

57  {/*{{{*/
58 
59  int num_proc = IssmComm::GetSize();
60 
61  /*Copy dataset*/
62  Vertices* output=new Vertices();
63  output->sorted = this->sorted;
64  output->numsorted = this->numsorted;
65  output->presorted = this->presorted;
66  for(vector<Object*>::iterator obj=this->objects.begin() ; obj < this->objects.end(); obj++ ) output->AddObject((*obj)->copy());
67 
68  /*Build id_offsets and sorted_ids*/
69  output->id_offsets=NULL;
70  output->sorted_ids=NULL;
71  int objsize = this->numsorted;
72  if(this->sorted && objsize>0 && this->id_offsets){
73  output->id_offsets=xNew<int>(objsize);
74  xMemCpy<int>(output->id_offsets,this->id_offsets,objsize);
75  }
76  if(this->sorted && objsize>0 && this->sorted_ids){
77  output->sorted_ids=xNew<int>(objsize);
78  xMemCpy<int>(output->sorted_ids,this->sorted_ids,objsize);
79  }
80 
81  /*Copy other fields*/
82  output->numberofvertices = this->numberofvertices;
85 
86  if(this->common_recv){
87  output->common_recv=xNew<int>(num_proc);
88  for(int i=0;i<num_proc;i++) output->common_recv[i]=this->common_recv[i];
89  }
90  if(this->common_send){
91  output->common_send=xNew<int>(num_proc);
92  for(int i=0;i<num_proc;i++) output->common_send[i]=this->common_send[i];
93  }
94  if(this->common_recv_ids){
95  output->common_recv_ids = xNew<int*>(num_proc);
96  for(int i=0;i<num_proc;i++){
97  output->common_recv_ids[i]=xNew<int>(this->common_recv[i]);
98  for(int j=0;j<this->common_recv[i];j++) output->common_recv_ids[i][j]=this->common_recv_ids[i][j];
99  }
100  }
101  if(this->common_send_ids){
102  output->common_send_ids = xNew<int*>(num_proc);
103  for(int i=0;i<num_proc;i++){
104  output->common_send_ids[i]=xNew<int>(this->common_send[i]);
105  for(int j=0;j<this->common_send[i];j++) output->common_send_ids[i][j]=this->common_send_ids[i][j];
106  }
107  }
108 
109  return output;
110 }

◆ Marshall()

void Vertices::Marshall ( char **  pmarshalled_data,
int *  pmarshalled_data_size,
int  marshall_direction 
)

Definition at line 112 of file Vertices.cpp.

112  { /*{{{*/
113 
114  int num_procs=IssmComm::GetSize();
115  int test = num_procs;
120  MARSHALLING(test);
121  if(test!=num_procs) _error_("number of cores is not the same as before");
122  MARSHALLING_DYNAMIC(this->common_recv,int,num_procs);
123  MARSHALLING_DYNAMIC(this->common_send,int,num_procs);
124  if(marshall_direction == MARSHALLING_BACKWARD){
125  this->common_recv_ids = xNew<int*>(num_procs);
126  this->common_send_ids = xNew<int*>(num_procs);
127  }
128  for(int i=0;i<num_procs;i++){
129  MARSHALLING_DYNAMIC(this->common_recv_ids[i],int,this->common_recv[i]);
130  MARSHALLING_DYNAMIC(this->common_send_ids[i],int,this->common_send[i]);
131  }
132  DataSet::Marshall(pmarshalled_data,pmarshalled_data_size,marshall_direction);
133 }

◆ Finalize()

void Vertices::Finalize ( IoModel iomodel)

Definition at line 173 of file Vertices.cpp.

173  {/*{{{*/
174 
175  /*Here we do 3 things:
176  * - count all vertices once for all so that we do not need to call MPI
177  * every time we need to know the total number of vertices
178  * - Distribute lids (local ids): masters first, slaves second
179  * - Distribute pids (parallel ids)
180  * */
181 
182  /*recover my_rank:*/
183  ISSM_MPI_Status status;
184  int my_rank = IssmComm::GetRank();
185  int num_procs = IssmComm::GetSize();
186 
187  /*1. set number of vertices once for all*/
188  this->numberofvertices_local=this->Size();
189  this->numberofmasters_local=0;
190  for(int i=0;i<this->Size();i++){
191  Vertex* vertex=xDynamicCast<Vertex*>(this->GetObjectByOffset(i));
192  if(!vertex->clone) this->numberofmasters_local++;
193  }
195 
196  /*2. Distribute lids (First: masters, then clones)*/
197  iomodel->my_vertices_lids=xNew<int>(this->numberofvertices);
198  for(int i=0;i<this->numberofvertices;i++) iomodel->my_vertices_lids[i] = -1;
199 
200  int lid = 0;
201  for(int i=0;i<this->Size();i++){
202  Vertex* vertex=xDynamicCast<Vertex*>(this->GetObjectByOffset(i));
203  if(!vertex->clone){
204  vertex->lid=lid;
205  iomodel->my_vertices_lids[vertex->sid] = lid;
206  lid++;
207  }
208  }
209  for(int i=0;i<this->Size();i++){
210  Vertex* vertex=xDynamicCast<Vertex*>(this->GetObjectByOffset(i));
211  if(vertex->clone){
212  vertex->lid=lid;
213  iomodel->my_vertices_lids[vertex->sid] = lid;
214  lid++;
215  }
216  }
217 
218  /*3. Distribute pids based on lids and offsets*/
219  int* all_num_masters=xNew<int>(num_procs);
220  ISSM_MPI_Gather(&this->numberofmasters_local,1,ISSM_MPI_INT,all_num_masters,1,ISSM_MPI_INT,0,IssmComm::GetComm());
221  ISSM_MPI_Bcast(all_num_masters,num_procs,ISSM_MPI_INT,0,IssmComm::GetComm());
222  int offset=0;
223  for(int i=0;i<my_rank;i++) offset+=all_num_masters[i];
224  xDelete<int>(all_num_masters);
225 
226  for(int i=0;i<this->Size();i++){
227  Vertex* vertex=xDynamicCast<Vertex*>(this->GetObjectByOffset(i));
228  vertex->pid = vertex->lid+offset;
229  }
230 
231  /* Share pids of masters and update pids of clones*/
232  int* truepids = xNew<int>(this->Size()); //only one alloc
233  for(int rank=0;rank<num_procs;rank++){
234  if(this->common_send[rank]){
235  int numids = this->common_send[rank];
236  for(int i=0;i<numids;i++){
237  Vertex* vertex=xDynamicCast<Vertex*>(this->GetObjectByOffset(this->common_send_ids[rank][i]));
238  truepids[i] = vertex->pid;
239  }
240  ISSM_MPI_Send(truepids,numids,ISSM_MPI_INT,rank,0,IssmComm::GetComm());
241  }
242  }
243  for(int rank=0;rank<num_procs;rank++){
244  if(this->common_recv[rank]){
245  int numids = this->common_recv[rank];
246  ISSM_MPI_Recv(truepids,numids,ISSM_MPI_INT,rank,0,IssmComm::GetComm(),&status);
247  for(int i=0;i<numids;i++){
248  Vertex* vertex=xDynamicCast<Vertex*>(this->GetObjectByOffset(this->common_recv_ids[rank][i]));
249  vertex->pid = truepids[i];
250  }
251  }
252  }
253  xDelete<int>(truepids);
254 }/*}}}*/

◆ NumberOfVertices()

int Vertices::NumberOfVertices ( void  )

Definition at line 255 of file Vertices.cpp.

255  {/*{{{*/
256  return this->numberofvertices;
257 }/*}}}*/

◆ NumberOfVerticesLocal()

int Vertices::NumberOfVerticesLocal ( void  )

Definition at line 258 of file Vertices.cpp.

258  {/*{{{*/
259  return this->numberofmasters_local;
260 }/*}}}*/

◆ NumberOfVerticesLocalAll()

int Vertices::NumberOfVerticesLocalAll ( void  )

Definition at line 261 of file Vertices.cpp.

261  {/*{{{*/
262  return this->numberofvertices_local;
263 }/*}}}*/

◆ LatLonList()

void Vertices::LatLonList ( IssmDouble **  lat,
IssmDouble **  lon 
)

Definition at line 135 of file Vertices.cpp.

135  {/*{{{*/
136 
137  /*output: */
138  IssmDouble* xyz_serial=NULL;
139 
140  /*recover my_rank:*/
141  int my_rank=IssmComm::GetRank();
142 
143  /*First, figure out number of vertices: */
144  int num_vertices=this->NumberOfVertices();
145 
146  /*Now, allocate vectors*/
147  Vector<IssmDouble>* lat = new Vector<IssmDouble>(num_vertices);
148  Vector<IssmDouble>* lon = new Vector<IssmDouble>(num_vertices);
149 
150  /*Go through vertices, and for each vertex, object, report it cpu: */
151  for(int i=0;i<this->Size();i++){
152  Vertex* vertex=xDynamicCast<Vertex*>(this->GetObjectByOffset(i));
153  lat->SetValue(vertex->sid,vertex->GetLatitude() ,INS_VAL);
154  lon->SetValue(vertex->sid,vertex->GetLongitude(),INS_VAL);
155  }
156 
157  /*Assemble:*/
158  lat->Assemble();
159  lon->Assemble();
160 
161  /*gather on cpu 0: */
162  IssmDouble* lat_serial=lat->ToMPISerial();
163  IssmDouble* lon_serial=lon->ToMPISerial();
164 
165  /*free ressources: */
166  *plat = lat_serial;
167  *plon = lon_serial;
168  delete lat;
169  delete lon;
170 }

Field Documentation

◆ numberofvertices

int Vertices::numberofvertices
private

Definition at line 18 of file Vertices.h.

◆ numberofvertices_local

int Vertices::numberofvertices_local
private

Definition at line 19 of file Vertices.h.

◆ numberofmasters_local

int Vertices::numberofmasters_local
private

Definition at line 20 of file Vertices.h.

◆ common_recv

int* Vertices::common_recv

Definition at line 22 of file Vertices.h.

◆ common_recv_ids

int** Vertices::common_recv_ids

Definition at line 23 of file Vertices.h.

◆ common_send

int* Vertices::common_send

Definition at line 24 of file Vertices.h.

◆ common_send_ids

int** Vertices::common_send_ids

Definition at line 25 of file Vertices.h.


The documentation for this class was generated from the following files:
DataSet::Size
int Size()
Definition: DataSet.cpp:399
Vertices::common_recv
int * common_recv
Definition: Vertices.h:22
Vertices
Declaration of Vertices class.
Definition: Vertices.h:15
Vertex::GetLongitude
IssmDouble GetLongitude(void)
Definition: Vertex.cpp:144
IssmDouble
double IssmDouble
Definition: types.h:37
Vertices::numberofvertices_local
int numberofvertices_local
Definition: Vertices.h:19
Vertices::common_send
int * common_send
Definition: Vertices.h:24
ISSM_MPI_Allreduce
int ISSM_MPI_Allreduce(void *sendbuf, void *recvbuf, int count, ISSM_MPI_Datatype datatype, ISSM_MPI_Op op, ISSM_MPI_Comm comm)
Definition: issmmpi.cpp:116
DataSet::AddObject
int AddObject(Object *object)
Definition: DataSet.cpp:252
Vertices::numberofvertices
int numberofvertices
Definition: Vertices.h:18
ISSM_MPI_SUM
#define ISSM_MPI_SUM
Definition: issmmpi.h:134
Vertex::clone
bool clone
Definition: Vertex.h:22
MARSHALLING_ENUM
#define MARSHALLING_ENUM(EN)
Definition: Marshalling.h:14
IssmComm::GetComm
static ISSM_MPI_Comm GetComm(void)
Definition: IssmComm.cpp:30
IoModel::my_vertices_lids
int * my_vertices_lids
Definition: IoModel.h:73
DataSet::objects
std::vector< Object * > objects
Definition: DataSet.h:19
Vertices::NumberOfVertices
int NumberOfVertices(void)
Definition: Vertices.cpp:255
DataSet::sorted_ids
int * sorted_ids
Definition: DataSet.h:28
Vertices::common_recv_ids
int ** common_recv_ids
Definition: Vertices.h:23
MARSHALLING_DYNAMIC
#define MARSHALLING_DYNAMIC(FIELD, TYPE, SIZE)
Definition: Marshalling.h:61
Vertices::Vertices
Vertices()
Definition: Vertices.cpp:26
Vertex::pid
int pid
Definition: Vertex.h:26
ISSM_MPI_INT
#define ISSM_MPI_INT
Definition: issmmpi.h:127
DataSet::sorted
int sorted
Definition: DataSet.h:25
IssmComm::GetSize
static int GetSize(void)
Definition: IssmComm.cpp:46
VerticesEnum
@ VerticesEnum
Definition: EnumDefinitions.h:1326
MARSHALLING
#define MARSHALLING(FIELD)
Definition: Marshalling.h:29
ISSM_MPI_Send
int ISSM_MPI_Send(void *buf, int count, ISSM_MPI_Datatype datatype, int dest, int tag, ISSM_MPI_Comm comm)
Definition: issmmpi.cpp:484
INS_VAL
@ INS_VAL
Definition: toolkitsenums.h:14
ISSM_MPI_Status
int ISSM_MPI_Status
Definition: issmmpi.h:121
DataSet::enum_type
int enum_type
Definition: DataSet.h:22
ISSM_MPI_Bcast
int ISSM_MPI_Bcast(void *buffer, int count, ISSM_MPI_Datatype datatype, int root, ISSM_MPI_Comm comm)
Definition: issmmpi.cpp:162
Vector::Assemble
void Assemble(void)
Definition: Vector.h:142
MARSHALLING_BACKWARD
@ MARSHALLING_BACKWARD
Definition: Marshalling.h:10
IssmComm::GetRank
static int GetRank(void)
Definition: IssmComm.cpp:34
_error_
#define _error_(StreamArgs)
Definition: exceptions.h:49
DataSet::GetObjectByOffset
Object * GetObjectByOffset(int offset)
Definition: DataSet.cpp:334
Vertex::GetLatitude
IssmDouble GetLatitude(void)
Definition: Vertex.cpp:140
Vertex::lid
int lid
Definition: Vertex.h:27
DataSet::numsorted
int numsorted
Definition: DataSet.h:27
Vertex
Definition: Vertex.h:19
ISSM_MPI_Recv
int ISSM_MPI_Recv(void *buf, int count, ISSM_MPI_Datatype datatype, int source, int tag, ISSM_MPI_Comm comm, ISSM_MPI_Status *status)
Definition: issmmpi.cpp:342
Vector::ToMPISerial
doubletype * ToMPISerial(void)
Definition: Vector.h:277
DataSet::id_offsets
int * id_offsets
Definition: DataSet.h:29
Vertices::numberofmasters_local
int numberofmasters_local
Definition: Vertices.h:20
Vector< IssmDouble >
Vector::SetValue
void SetValue(int dof, doubletype value, InsMode mode)
Definition: Vector.h:163
Vertex::sid
int sid
Definition: Vertex.h:25
ISSM_MPI_Gather
int ISSM_MPI_Gather(void *sendbuf, int sendcnt, ISSM_MPI_Datatype sendtype, void *recvbuf, int recvcnt, ISSM_MPI_Datatype recvtype, int root, ISSM_MPI_Comm comm)
Definition: issmmpi.cpp:242
Vertices::common_send_ids
int ** common_send_ids
Definition: Vertices.h:25
DataSet::presorted
int presorted
Definition: DataSet.h:26
DataSet::Marshall
void Marshall(char **pmarshalled_data, int *pmarshalled_data_size, int marshall_direction)
Definition: DataSet.cpp:93