 |
Ice Sheet System Model
4.18
Code documentation
|
Declaration of Nodes class.
More...
#include <Nodes.h>
Declaration of Nodes class.
Declaration of Nodes class. Nodes are vector lists of objects (Containers) of Node objects. Node objects are the degrees of freedom (DOFs) for a particular analysis type (not to be confused with a vertex, which defines the (x,y,z) location of a point).
Definition at line 19 of file Nodes.h.
◆ Nodes()
◆ ~Nodes()
Definition at line 30 of file Nodes.cpp.
36 for(
int i=0;i<num_proc;i++) if(common_recv_ids[i]) xDelete<int>(
common_recv_ids[i]);
40 for(
int i=0;i<num_proc;i++) if(common_send_ids[i]) xDelete<int>(
common_send_ids[i]);
◆ Copy()
Nodes * Nodes::Copy |
( |
void |
| ) |
|
Definition at line 48 of file Nodes.cpp.
57 for(vector<Object*>::iterator obj=this->
objects.begin() ; obj < this->
objects.end(); obj++ ){
67 xMemCpy<int>(output->
id_offsets,this->id_offsets,objsize);
71 xMemCpy<int>(output->
sorted_ids,this->sorted_ids,objsize);
81 for(
int i=0;i<num_proc;i++) output->
common_recv[i]=this->common_recv[i];
85 for(
int i=0;i<num_proc;i++) output->
common_send[i]=this->common_send[i];
89 for(
int i=0;i<num_proc;i++){
96 for(
int i=0;i<num_proc;i++){
◆ Marshall()
void Nodes::Marshall |
( |
char ** |
pmarshalled_data, |
|
|
int * |
pmarshalled_data_size, |
|
|
int |
marshall_direction |
|
) |
| |
Definition at line 105 of file Nodes.cpp.
108 int test = num_procs;
115 if(test!=num_procs)
_error_(
"number of cores is not the same as before");
122 for(
int i=0;i<num_procs;i++){
129 if(this->
Size()==0)
return;
133 for(
int i=0;i<num_procs;i++){
◆ DistributeDofs()
void Nodes::DistributeDofs |
( |
int |
SETENUM | ) |
|
Definition at line 139 of file Nodes.cpp.
147 for(
int i=0;i<this->
Size();i++){
154 for(
int i=0;i<this->
Size();i++){
159 int dofcount_local = dofcount;
160 for(
int i=0;i<this->
Size();i++){
169 int* alldofcount=xNew<int>(num_procs);
175 for(
int i=0;i<my_rank;i++) offset+=alldofcount[i];
176 xDelete<int>(alldofcount);
178 for(
int i=0;i<this->
Size();i++){
186 int maxdofspernode = this->
MaxNumDofs(setenum);
187 int* truedofs = xNewZeroInit<int>(this->
Size()*maxdofspernode);
188 for(
int rank=0;rank<num_procs;rank++){
191 for(
int i=0;i<numids;i++){
198 for(
int rank=0;rank<num_procs;rank++){
202 for(
int i=0;i<numids;i++){
208 xDelete<int>(truedofs);
211 for(
int i=0;i<this->
Size();i++){
◆ Finalize()
void Nodes::Finalize |
( |
void |
| ) |
|
Definition at line 217 of file Nodes.cpp.
235 for(
int i=0;i<this->
Size();i++){
243 for(
int i=0;i<this->
Size();i++){
247 for(
int i=0;i<this->
Size();i++){
253 int* all_num_masters=xNew<int>(num_procs);
257 for(
int i=0;i<my_rank;i++) offset+=all_num_masters[i];
258 xDelete<int>(all_num_masters);
260 for(
int i=0;i<this->
Size();i++){
262 node->
pid = node->
lid+offset;
266 int* truepids = xNew<int>(this->
Size());
267 for(
int rank=0;rank<num_procs;rank++){
270 for(
int i=0;i<numids;i++){
272 truepids[i] = node->
pid;
277 for(
int rank=0;rank<num_procs;rank++){
281 for(
int i=0;i<numids;i++){
283 node->
pid = truepids[i];
287 xDelete<int>(truepids);
◆ MaxNumDofs()
int Nodes::MaxNumDofs |
( |
int |
setenum | ) |
|
Definition at line 294 of file Nodes.cpp.
300 for(
int i=0;i<this->
Size();i++){
304 if(numdofs>
max)
max=numdofs;
◆ NumberOfDofs()
int Nodes::NumberOfDofs |
( |
int |
setenum | ) |
|
◆ NumberOfDofsLocal()
int Nodes::NumberOfDofsLocal |
( |
int |
setenum | ) |
|
◆ NumberOfDofsLocalAll()
int Nodes::NumberOfDofsLocalAll |
( |
int |
setenum | ) |
|
◆ NumberOfNodes()
int Nodes::NumberOfNodes |
( |
void |
| ) |
|
◆ NumberOfNodesLocal()
int Nodes::NumberOfNodesLocal |
( |
void |
| ) |
|
◆ NumberOfNodesLocalAll()
int Nodes::NumberOfNodesLocalAll |
( |
void |
| ) |
|
◆ RequiresDofReindexing()
bool Nodes::RequiresDofReindexing |
( |
void |
| ) |
|
Definition at line 369 of file Nodes.cpp.
375 for(
int i=0;i<this->
Size();i++){
386 if(allflag)
return true;
◆ CheckDofListAcrossPartitions()
void Nodes::CheckDofListAcrossPartitions |
( |
void |
| ) |
|
Definition at line 391 of file Nodes.cpp.
407 for(
int i=0;i<this->
Size();i++){
415 for(
int j=0;j<node->
gsize;j++){
417 if(node->
s_set[j])
_error_(
"a degree of freedom is both in f and s set!");
422 if(node->
s_set[j]==0)
_error_(
"a degree of freedom is neither in f nor in s set!");
435 for(
int i=0;i<this->
Size();i++){
442 for(
int j=0;j<node->
gsize;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]);
451 if(local_dofs_check[index] != -1.){
452 _error_(
"Dof #"<<j<<
" of node sid "<<node->
Sid()<<
" not consistently in s set");
461 xDelete<IssmDouble>(local_dofs_check);
◆ GetLocalVectorWithClonesGset()
Definition at line 463 of file Nodes.cpp.
476 int *indices_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);
483 IssmDouble *local_ug = xNew<IssmDouble>(glocalsize);
484 xMemCpy<IssmDouble>(local_ug,local_ug_masters,glocalsize_masters);
485 xDelete<IssmDouble>(local_ug_masters);
489 IssmDouble* buffer = xNew<IssmDouble>(this->
Size()*maxdofspernode,
"t");
491 IssmDouble* buffer = xNew<IssmDouble>(this->
Size()*maxdofspernode);
493 for(
int rank=0;rank<num_procs;rank++){
496 for(
int i=0;i<numids;i++){
500 for(
int j=0;j<node->
gsize;j++) buffer[i*maxdofspernode+j]=local_ug[node->
gdoflist_local[j]];
505 for(
int rank=0;rank<num_procs;rank++){
509 for(
int i=0;i<numids;i++){
512 for(
int j=0;j<node->
gsize;j++) local_ug[node->
gdoflist_local[j]] = buffer[i*maxdofspernode+j];
516 xDelete<IssmDouble>(buffer);
519 *plocal_ug = local_ug;
◆ numberofnodes
◆ numberofnodes_local
int Nodes::numberofnodes_local |
|
private |
◆ numberofmasters_local
int Nodes::numberofmasters_local |
|
private |
◆ common_recv
◆ common_recv_ids
int** Nodes::common_recv_ids |
◆ common_send
◆ common_send_ids
int** Nodes::common_send_ids |
The documentation for this class was generated from the following files:
int numberofmasters_local
void UpdateCloneDofs(int *alltruerows, int setenum)
Declaration of Nodes class.
#define _printf0_(StreamArgs)
int ISSM_MPI_Allreduce(void *sendbuf, void *recvbuf, int count, ISSM_MPI_Datatype datatype, ISSM_MPI_Op op, ISSM_MPI_Comm comm)
int AddObject(Object *object)
int MaxNumDofs(int setenum)
#define MARSHALLING_ENUM(EN)
static ISSM_MPI_Comm GetComm(void)
std::vector< Object * > objects
void DistributeLocalDofs(int *pdofcount, int setenum)
void DistributeGlobalDofsMasters(int dofcount, int setenum)
#define MARSHALLING_DYNAMIC(FIELD, TYPE, SIZE)
int NumberOfDofsLocal(int setenum)
void ReindexingDone(void)
void GetLocalVectorWithClonesGset(IssmDouble **plocal_vector, Vector< IssmDouble > *vector)
#define MARSHALLING(FIELD)
int ISSM_MPI_Send(void *buf, int count, ISSM_MPI_Datatype datatype, int dest, int tag, ISSM_MPI_Comm comm)
void AllocateDofLists(int setenum)
int ISSM_MPI_Bcast(void *buffer, int count, ISSM_MPI_Datatype datatype, int root, ISSM_MPI_Comm comm)
int NumberOfDofsLocalAll(int setenum)
#define _error_(StreamArgs)
Object * GetObjectByOffset(int offset)
int ISSM_MPI_Recv(void *buf, int count, ISSM_MPI_Datatype datatype, int source, int tag, ISSM_MPI_Comm comm, ISSM_MPI_Status *status)
IssmDouble max(IssmDouble a, IssmDouble b)
int NumberOfDofs(int setenum)
int GetNumberOfDofs(int approximation_enum, int setenum)
bool RequiresDofReindexing(void)
void SetValue(int dof, doubletype value, InsMode mode)
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)
void ShowMasterDofs(int *truerows, int setenum)
void Marshall(char **pmarshalled_data, int *pmarshalled_data_size, int marshall_direction)