Ice Sheet System Model  4.18
Code documentation
Nodes.cpp
Go to the documentation of this file.
1 /*
2  * \file Nodes.cpp
3  * \brief: Implementation of Nodes class, derived from DataSet class.
4  */
5 
6 /*Headers*/
7 #ifdef HAVE_CONFIG_H
8  #include <config.h>
9 #else
10 #error "Cannot compile with HAVE_CONFIG_H symbol! run configure first!"
11 #endif
12 #include "../shared/io/Comm/IssmComm.h"
13 #include "./Nodes.h"
14 #include "./Node.h"
15 using namespace std;
16 
17 /*Object constructors and destructor*/
18 Nodes::Nodes(){/*{{{*/
19  this->enum_type = NodesEnum;
20  this->common_recv = NULL;
21  this->common_recv_ids = NULL;
22  this->common_send = NULL;
23  this->common_send_ids = NULL;
24  this->numberofnodes = -1;
25  this->numberofnodes_local = -1;
26  this->numberofmasters_local = -1;
27  return;
28 }
29 /*}}}*/
30 Nodes::~Nodes(){/*{{{*/
31  int num_proc=IssmComm::GetSize();
32 
33  if(this->common_recv) xDelete<int>(common_recv);
34  if(this->common_send) xDelete<int>(common_send);
35  if(this->common_recv_ids){
36  for(int i=0;i<num_proc;i++) if(common_recv_ids[i]) xDelete<int>(common_recv_ids[i]);
37  xDelete<int*>(common_recv_ids);
38  }
39  if(this->common_send_ids){
40  for(int i=0;i<num_proc;i++) if(common_send_ids[i]) xDelete<int>(common_send_ids[i]);
41  xDelete<int*>(common_send_ids);
42  }
43  return;
44 }
45 /*}}}*/
46 
47 /*Numerics*/
48 Nodes* Nodes::Copy() {/*{{{*/
49 
50  int num_proc = IssmComm::GetSize();
51 
52  /*Copy dataset*/
53  Nodes* output=new Nodes();
54  output->sorted=this->sorted;
55  output->numsorted=this->numsorted;
56  output->presorted=this->presorted;
57  for(vector<Object*>::iterator obj=this->objects.begin() ; obj < this->objects.end(); obj++ ){
58  output->AddObject((*obj)->copy());
59  }
60 
61  /*Build id_offsets and sorted_ids*/
62  int objsize = this->numsorted;
63  output->id_offsets=NULL;
64  output->sorted_ids=NULL;
65  if(this->sorted && objsize>0 && this->id_offsets){
66  output->id_offsets=xNew<int>(objsize);
67  xMemCpy<int>(output->id_offsets,this->id_offsets,objsize);
68  }
69  if(this->sorted && objsize>0 && this->sorted_ids){
70  output->sorted_ids=xNew<int>(objsize);
71  xMemCpy<int>(output->sorted_ids,this->sorted_ids,objsize);
72  }
73 
74  /*Copy other fields*/
75  output->numberofnodes = this->numberofnodes;
76  output->numberofnodes_local = this->numberofnodes_local;
77  output->numberofmasters_local = this->numberofmasters_local;
78 
79  if(this->common_recv){
80  output->common_recv=xNew<int>(num_proc);
81  for(int i=0;i<num_proc;i++) output->common_recv[i]=this->common_recv[i];
82  }
83  if(this->common_send){
84  output->common_send=xNew<int>(num_proc);
85  for(int i=0;i<num_proc;i++) output->common_send[i]=this->common_send[i];
86  }
87  if(this->common_recv_ids){
88  output->common_recv_ids = xNew<int*>(num_proc);
89  for(int i=0;i<num_proc;i++){
90  output->common_recv_ids[i]=xNew<int>(this->common_recv[i]);
91  for(int j=0;j<this->common_recv[i];j++) output->common_recv_ids[i][j]=this->common_recv_ids[i][j];
92  }
93  }
94  if(this->common_send_ids){
95  output->common_send_ids = xNew<int*>(num_proc);
96  for(int i=0;i<num_proc;i++){
97  output->common_send_ids[i]=xNew<int>(this->common_send[i]);
98  for(int j=0;j<this->common_send[i];j++) output->common_send_ids[i][j]=this->common_send_ids[i][j];
99  }
100  }
101 
102  return output;
103 }
104 /*}}}*/
105 void Nodes::Marshall(char** pmarshalled_data,int* pmarshalled_data_size, int marshall_direction){ /*{{{*/
106 
107  int num_procs=IssmComm::GetSize();
108  int test = num_procs;
110  MARSHALLING(numberofnodes);
111  MARSHALLING(numberofnodes_local);
112  MARSHALLING(numberofmasters_local);
113 
114  MARSHALLING(test);
115  if(test!=num_procs) _error_("number of cores is not the same as before");
116 
117  DataSet::Marshall(pmarshalled_data,pmarshalled_data_size,marshall_direction);
118 
119  if(marshall_direction == MARSHALLING_BACKWARD){
120  this->common_recv_ids = xNew<int*>(num_procs);
121  this->common_send_ids = xNew<int*>(num_procs);
122  for(int i=0;i<num_procs;i++){
123  this->common_recv_ids[i] = NULL;
124  this->common_send_ids[i] = NULL;
125  }
126  }
127 
128  /*Stop here if no nodes*/
129  if(this->Size()==0) return;
130 
131  MARSHALLING_DYNAMIC(this->common_recv,int,num_procs);
132  MARSHALLING_DYNAMIC(this->common_send,int,num_procs);
133  for(int i=0;i<num_procs;i++){
134  if(this->common_recv[i]) MARSHALLING_DYNAMIC(this->common_recv_ids[i],int,this->common_recv[i]);
135  if(this->common_send[i]) MARSHALLING_DYNAMIC(this->common_send_ids[i],int,this->common_send[i]);
136  }
137 }
138 /*}}}*/
139 void Nodes::DistributeDofs(int setenum){/*{{{*/
140 
141  /*recover my_rank:*/
142  ISSM_MPI_Status status;
143  int my_rank = IssmComm::GetRank();
144  int num_procs = IssmComm::GetSize();
145 
146  /*First, allocate dof lists based on fset and sset*/
147  for(int i=0;i<this->Size();i++){
148  Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i));
149  node->AllocateDofLists(setenum);
150  }
151 
152  /*Now, Build local dofs for masters first*/
153  int dofcount=0;
154  for(int i=0;i<this->Size();i++){
155  Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i));
156  if(!node->IsClone()) node->DistributeLocalDofs(&dofcount,setenum);
157  }
158  /*Build local dofs for clones, they always will be at the end*/
159  int dofcount_local = dofcount;
160  for(int i=0;i<this->Size();i++){
161  Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i));
162  if(node->IsClone()) node->DistributeLocalDofs(&dofcount_local,setenum);
163  }
164 
165  /* Now every object has distributed dofs, but locally, and with a dof count starting from
166  * 0. This means the dofs between all the cpus are not unique. We now offset the dofs of each
167  * cpus by the total last (master) dofs of the previus cpu, starting from 0.
168  * First: get number of dofs for each cpu*/
169  int* alldofcount=xNew<int>(num_procs);
170  ISSM_MPI_Gather(&dofcount,1,ISSM_MPI_INT,alldofcount,1,ISSM_MPI_INT,0,IssmComm::GetComm());
171  ISSM_MPI_Bcast(alldofcount,num_procs,ISSM_MPI_INT,0,IssmComm::GetComm());
172 
173  /* Every cpu should start its own dof count at the end of the dofcount from cpu-1*/
174  int offset=0;
175  for(int i=0;i<my_rank;i++) offset+=alldofcount[i];
176  xDelete<int>(alldofcount);
177 
178  for(int i=0;i<this->Size();i++){
179  Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i));
180  node->DistributeGlobalDofsMasters(offset,setenum);
181  }
182 
183  /* Finally, remember that cpus may have skipped some objects, because they were clones. For every
184  * object that is not a clone, tell them to show their dofs, so that later on, they can get picked
185  * up by their clones: */
186  int maxdofspernode = this->MaxNumDofs(setenum);
187  int* truedofs = xNewZeroInit<int>(this->Size()*maxdofspernode); //only one alloc
188  for(int rank=0;rank<num_procs;rank++){
189  if(this->common_send[rank]){
190  int numids = this->common_send[rank];
191  for(int i=0;i<numids;i++){
192  Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(this->common_send_ids[rank][i]));
193  node->ShowMasterDofs(&truedofs[i*maxdofspernode+0],setenum);
194  }
195  ISSM_MPI_Send(truedofs,numids*maxdofspernode,ISSM_MPI_INT,rank,0,IssmComm::GetComm());
196  }
197  }
198  for(int rank=0;rank<num_procs;rank++){
199  if(this->common_recv[rank]){
200  int numids = this->common_recv[rank];
201  ISSM_MPI_Recv(truedofs,numids*maxdofspernode,ISSM_MPI_INT,rank,0,IssmComm::GetComm(),&status);
202  for(int i=0;i<numids;i++){
203  Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(this->common_recv_ids[rank][i]));
204  node->UpdateCloneDofs(&truedofs[i*maxdofspernode+0],setenum);
205  }
206  }
207  }
208  xDelete<int>(truedofs);
209 
210  /*Update indexingupdateflag*/
211  for(int i=0;i<this->Size();i++){
212  Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i));
213  node->ReindexingDone();
214  }
215 }
216 /*}}}*/
217 void Nodes::Finalize(){/*{{{*/
218 
219  /*Here we do 4 things:
220  * - count all nodes once for all so that we do not need to call MPI
221  * every time we need to know the total number of vertices
222  * - Distribute lids (local ids): masters first, slaves second
223  * - Distribute pids (parallel ids)
224  * - Distribute Gset once for all
225  */
226 
227  /*recover my_rank:*/
228  ISSM_MPI_Status status;
229  int my_rank = IssmComm::GetRank();
230  int num_procs = IssmComm::GetSize();
231 
232  /*1. set number of nodes once for all*/
233  this->numberofnodes_local=this->Size();
234  this->numberofmasters_local=0;
235  for(int i=0;i<this->Size();i++){
236  Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i));
237  if(!node->clone) this->numberofmasters_local++;
238  }
239  ISSM_MPI_Allreduce((void*)&this->numberofmasters_local,(void*)&this->numberofnodes,1,ISSM_MPI_INT,ISSM_MPI_SUM,IssmComm::GetComm());
240 
241  /*2. Distribute lids (First: masters, then clones)*/
242  int lid = 0;
243  for(int i=0;i<this->Size();i++){
244  Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i));
245  if(!node->clone) node->lid=lid++;
246  }
247  for(int i=0;i<this->Size();i++){
248  Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i));
249  if(node->clone) node->lid=lid++;
250  }
251 
252  /*3. Distribute pids based on lids and offsets*/
253  int* all_num_masters=xNew<int>(num_procs);
254  ISSM_MPI_Gather(&this->numberofmasters_local,1,ISSM_MPI_INT,all_num_masters,1,ISSM_MPI_INT,0,IssmComm::GetComm());
255  ISSM_MPI_Bcast(all_num_masters,num_procs,ISSM_MPI_INT,0,IssmComm::GetComm());
256  int offset=0;
257  for(int i=0;i<my_rank;i++) offset+=all_num_masters[i];
258  xDelete<int>(all_num_masters);
259 
260  for(int i=0;i<this->Size();i++){
261  Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i));
262  node->pid = node->lid+offset;
263  }
264 
265  /* Share pids of masters and update pids of clones*/
266  int* truepids = xNew<int>(this->Size()); //only one alloc
267  for(int rank=0;rank<num_procs;rank++){
268  if(this->common_send[rank]){
269  int numids = this->common_send[rank];
270  for(int i=0;i<numids;i++){
271  Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(this->common_send_ids[rank][i]));
272  truepids[i] = node->pid;
273  }
274  ISSM_MPI_Send(truepids,numids,ISSM_MPI_INT,rank,0,IssmComm::GetComm());
275  }
276  }
277  for(int rank=0;rank<num_procs;rank++){
278  if(this->common_recv[rank]){
279  int numids = this->common_recv[rank];
280  ISSM_MPI_Recv(truepids,numids,ISSM_MPI_INT,rank,0,IssmComm::GetComm(),&status);
281  for(int i=0;i<numids;i++){
282  Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(this->common_recv_ids[rank][i]));
283  node->pid = truepids[i];
284  }
285  }
286  }
287  xDelete<int>(truepids);
288 
289  /*4. Distribute G dofs once for all*/
290  //this->DistributeDofs(GsetEnum);
291 
292  return;
293 }/*}}}*/
294 int Nodes::MaxNumDofs(int setenum){/*{{{*/
295 
296  int max=0;
297  int allmax;
298 
299  /*Now go through all nodes, and get how many dofs they own, unless they are clone nodes: */
300  for(int i=0;i<this->Size();i++){
301  Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i));
302 
303  int numdofs=node->GetNumberOfDofs(NoneApproximationEnum,setenum);
304  if(numdofs>max) max=numdofs;
305  }
306 
307  /*Grab max of all cpus: */
309  max=allmax;
310 
311  return max;
312 }
313 /*}}}*/
314 int Nodes::NumberOfDofs(int setenum){/*{{{*/
315 
316  int allnumdofs;
317 
318  /*Get number of dofs on current cpu (excluding clones)*/
319  int numdofs=this->NumberOfDofsLocal(setenum);
320 
321  /*Gather from all cpus: */
322  ISSM_MPI_Allreduce ((void*)&numdofs,(void*)&allnumdofs,1,ISSM_MPI_INT,ISSM_MPI_SUM,IssmComm::GetComm());
323  return allnumdofs;
324 }
325 /*}}}*/
326 int Nodes::NumberOfDofsLocal(int setenum){/*{{{*/
327 
328  int numdofs=0;
329 
330  /*Now go through all nodes, and get how many dofs they own, unless they are clone nodes: */
331  for(int i=0;i<this->Size();i++){
332  Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i));
333 
334  /*Ok, this object is a node, ask it to plug values into partition: */
335  if (!node->IsClone()){
336  numdofs+=node->GetNumberOfDofs(NoneApproximationEnum,setenum);
337  }
338  }
339 
340  return numdofs;
341 }
342 /*}}}*/
343 int Nodes::NumberOfDofsLocalAll(int setenum){/*{{{*/
344 
345  /*go through all nodes, and get how many dofs they own, unless they are clone nodes: */
346  int numdofs=0;
347  for(int i=0;i<this->Size();i++){
348  Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i));
349  numdofs+=node->GetNumberOfDofs(NoneApproximationEnum,setenum);
350  }
351  return numdofs;
352 }
353 /*}}}*/
354 int Nodes::NumberOfNodes(void){/*{{{*/
355 
356  return this->numberofnodes;
357 }
358 /*}}}*/
359 int Nodes::NumberOfNodesLocal(void){/*{{{*/
360 
361  return this->numberofmasters_local;
362 }
363 /*}}}*/
365 
366  return this->numberofnodes_local;
367 }
368 /*}}}*/
370 
371  int flag = 0;
372  int allflag;
373 
374  /*Now go through all nodes, and get how many dofs they own, unless they are clone nodes: */
375  for(int i=0;i<this->Size();i++){
376  Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i));
377  if(node->RequiresDofReindexing()){
378  flag = 1;
379  break;
380  }
381  }
382 
383  /*Grab max of all cpus: */
384  ISSM_MPI_Allreduce((void*)&flag,(void*)&allflag,1,ISSM_MPI_INT,ISSM_MPI_MAX,IssmComm::GetComm());
385 
386  if(allflag) return true;
387  else return false;
388 }
389 /*}}}*/
390 
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 }/*}}}*/
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 }/*}}}*/
Nodes::numberofmasters_local
int numberofmasters_local
Definition: Nodes.h:24
Node::UpdateCloneDofs
void UpdateCloneDofs(int *alltruerows, int setenum)
Definition: Node.cpp:1012
Nodes::Nodes
Nodes()
Definition: Nodes.cpp:18
_assert_
#define _assert_(ignore)
Definition: exceptions.h:37
IssmDouble
double IssmDouble
Definition: types.h:37
Nodes
Declaration of Nodes class.
Definition: Nodes.h:19
_printf0_
#define _printf0_(StreamArgs)
Definition: Print.h:29
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
FemModel::GetLocalVectorWithClonesGset
void GetLocalVectorWithClonesGset(IssmDouble **plocal_ug, Vector< IssmDouble > *ug)
Definition: FemModel.cpp:1337
Nodes::RequiresDofReindexing
bool RequiresDofReindexing(void)
Definition: Nodes.cpp:369
DataSet::AddObject
int AddObject(Object *object)
Definition: DataSet.cpp:252
Nodes::common_send
int * common_send
Definition: Nodes.h:28
Nodes::MaxNumDofs
int MaxNumDofs(int setenum)
Definition: Nodes.cpp:294
ISSM_MPI_SUM
#define ISSM_MPI_SUM
Definition: issmmpi.h:134
MARSHALLING_ENUM
#define MARSHALLING_ENUM(EN)
Definition: Marshalling.h:14
IssmComm::GetComm
static ISSM_MPI_Comm GetComm(void)
Definition: IssmComm.cpp:30
FemModel::Size
int Size(void)
Definition: FemModel.cpp:638
Nodes::numberofnodes_local
int numberofnodes_local
Definition: Nodes.h:23
Nodes::NumberOfNodesLocalAll
int NumberOfNodesLocalAll(void)
Definition: Nodes.cpp:364
Vector::GetLocalVector
void GetLocalVector(doubletype **pvector, int **pindices)
Definition: Vector.h:219
Node::DistributeLocalDofs
void DistributeLocalDofs(int *pdofcount, int setenum)
Definition: Node.cpp:943
ISSM_MPI_MAX
#define ISSM_MPI_MAX
Definition: issmmpi.h:131
Node::clone
bool clone
Definition: Node.h:27
DataSet::sorted_ids
int * sorted_ids
Definition: DataSet.h:28
Nodes::common_recv_ids
int ** common_recv_ids
Definition: Nodes.h:27
VerboseModule
bool VerboseModule(void)
Definition: Verbosity.cpp:23
Node::DistributeGlobalDofsMasters
void DistributeGlobalDofsMasters(int dofcount, int setenum)
Definition: Node.cpp:970
Nodes::Marshall
void Marshall(char **pmarshalled_data, int *pmarshalled_data_size, int marshall_direction)
Definition: Nodes.cpp:105
Nodes::~Nodes
~Nodes()
Definition: Nodes.cpp:30
MARSHALLING_DYNAMIC
#define MARSHALLING_DYNAMIC(FIELD, TYPE, SIZE)
Definition: Marshalling.h:61
Node::s_set
bool * s_set
Definition: Node.h:54
ISSM_MPI_DOUBLE
#define ISSM_MPI_DOUBLE
Definition: issmmpi.h:125
Nodes::NumberOfNodes
int NumberOfNodes(void)
Definition: Nodes.cpp:354
Nodes::NumberOfDofsLocal
int NumberOfDofsLocal(int setenum)
Definition: Nodes.cpp:326
ISSM_MPI_INT
#define ISSM_MPI_INT
Definition: issmmpi.h:127
Nodes::common_send_ids
int ** common_send_ids
Definition: Nodes.h:29
NoneApproximationEnum
@ NoneApproximationEnum
Definition: EnumDefinitions.h:1201
Nodes::CheckDofListAcrossPartitions
void CheckDofListAcrossPartitions(void)
Definition: Nodes.cpp:391
DataSet::sorted
int sorted
Definition: DataSet.h:25
Node::ReindexingDone
void ReindexingDone(void)
Definition: Node.cpp:807
Node::pid
int pid
Definition: Node.h:31
NodesEnum
@ NodesEnum
Definition: EnumDefinitions.h:275
IssmComm::GetSize
static int GetSize(void)
Definition: IssmComm.cpp:46
GsetEnum
@ GsetEnum
Definition: EnumDefinitions.h:1093
Nodes::GetLocalVectorWithClonesGset
void GetLocalVectorWithClonesGset(IssmDouble **plocal_vector, Vector< IssmDouble > *vector)
Definition: Nodes.cpp:463
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
Nodes::common_recv
int * common_recv
Definition: Nodes.h:26
INS_VAL
@ INS_VAL
Definition: toolkitsenums.h:14
Nodes::numberofnodes
int numberofnodes
Definition: Nodes.h:22
ISSM_MPI_Status
int ISSM_MPI_Status
Definition: issmmpi.h:121
Node.h
: header file for node object
Node::IsClone
int IsClone()
Definition: Node.cpp:801
Node::AllocateDofLists
void AllocateDofLists(int setenum)
Definition: Node.cpp:905
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
Nodes::DistributeDofs
void DistributeDofs(int SETENUM)
Definition: Nodes.cpp:139
MARSHALLING_BACKWARD
@ MARSHALLING_BACKWARD
Definition: Marshalling.h:10
Nodes::Copy
Nodes * Copy()
Definition: Nodes.cpp:48
Nodes::Finalize
void Finalize(void)
Definition: Nodes.cpp:217
IssmComm::GetRank
static int GetRank(void)
Definition: IssmComm.cpp:34
Nodes::NumberOfDofsLocalAll
int NumberOfDofsLocalAll(int setenum)
Definition: Nodes.cpp:343
Node
Definition: Node.h:23
Node::Sid
int Sid(void)
Definition: Node.cpp:622
_error_
#define _error_(StreamArgs)
Definition: exceptions.h:49
Nodes.h
Node::f_set
bool * f_set
Definition: Node.h:53
Node::lid
int lid
Definition: Node.h:30
Node::gsize
int gsize
Definition: Node.h:44
Node::gdoflist_local
int * gdoflist_local
Definition: Node.h:64
DataSet::numsorted
int numsorted
Definition: DataSet.h:27
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
max
IssmDouble max(IssmDouble a, IssmDouble b)
Definition: extrema.cpp:24
Node::fdoflist
int * fdoflist
Definition: Node.h:62
DataSet::id_offsets
int * id_offsets
Definition: DataSet.h:29
Node::gdoflist
int * gdoflist
Definition: Node.h:61
Nodes::NumberOfDofs
int NumberOfDofs(int setenum)
Definition: Nodes.cpp:314
Vector< IssmDouble >
Node::GetNumberOfDofs
int GetNumberOfDofs(int approximation_enum, int setenum)
Definition: Node.cpp:741
Nodes::NumberOfNodesLocal
int NumberOfNodesLocal(void)
Definition: Nodes.cpp:359
Node::RequiresDofReindexing
bool RequiresDofReindexing(void)
Definition: Node.cpp:822
Vector::SetValue
void SetValue(int dof, doubletype value, InsMode mode)
Definition: Vector.h:163
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
DataSet::presorted
int presorted
Definition: DataSet.h:26
Node::ShowMasterDofs
void ShowMasterDofs(int *truerows, int setenum)
Definition: Node.cpp:991
DataSet::Marshall
void Marshall(char **pmarshalled_data, int *pmarshalled_data_size, int marshall_direction)
Definition: DataSet.cpp:93