Changeset 23627


Ignore:
Timestamp:
01/11/19 21:10:55 (6 years ago)
Author:
Mathieu Morlighem
Message:

CHG: revisited how dofs are allocated and distributed, now locally as well

Location:
issm/trunk-jpl/src/c/classes
Files:
3 edited

Legend:

Unmodified
Added
Removed
  • issm/trunk-jpl/src/c/classes/Node.cpp

    r23614 r23627  
    1818/*Node constructors and destructors:*/
    1919Node::Node(){/*{{{*/
    20         this->approximation=0;
    21 
    22         this->gsize    = -1;
    23         this->fsize    = -1;
    24         this->ssize    = -1;
    25         this->clone    = false;
    26         this->active   = true;
    27         this->freeze   = false;
    28         this->f_set    = NULL;
    29         this->s_set    = NULL;
    30         this->svalues  = NULL;
    31         this->doftype  = NULL;
    32         this->gdoflist = NULL;
    33         this->fdoflist = NULL;
    34         this->sdoflist = NULL;
     20        this->approximation  = 0;
     21        this->gsize          = -1;
     22        this->fsize          = -1;
     23        this->ssize          = -1;
     24        this->clone          = false;
     25        this->active         = true;
     26        this->freeze         = false;
     27        this->f_set          = NULL;
     28        this->s_set          = NULL;
     29        this->svalues        = NULL;
     30        this->doftype        = NULL;
     31        this->gdoflist       = NULL;
     32        this->fdoflist       = NULL;
     33        this->sdoflist       = NULL;
     34        this->gdoflist_local = NULL;
     35        this->fdoflist_local = NULL;
     36        this->sdoflist_local = NULL;
    3537}
    3638/*}}}*/
     
    4446        this->analysis_enum = node_analysis;
    4547        this->clone         = node_clone;
     48        this->active        = true;
     49        this->freeze        = false;
    4650
    4751        /*Initialize coord_system: Identity matrix by default*/
     
    4953        for(int k=0;k<3;k++) this->coord_system[k][k]=1.0;
    5054
    51         this->gsize    = -1;
    52         this->fsize    = -1;
    53         this->ssize    = -1;
    54         this->active   = true;
    55         this->freeze   = false;
    56         this->f_set    = NULL;
    57         this->s_set    = NULL;
    58         this->svalues  = NULL;
    59         this->doftype  = NULL;
    60         this->gdoflist = NULL;
    61         this->fdoflist = NULL;
    62         this->sdoflist = NULL;
     55        this->approximation=0;
     56        if(analysis_enum==StressbalanceAnalysisEnum) this->approximation=in_approximation;
    6357
    6458        /*indexing:*/
    6559        this->indexingupdate = true;
     60        this->doftype        = NULL;
    6661        Analysis *analysis = EnumToAnalysis(analysis_enum);
    67         int      *doftypes = NULL;
    6862        this->gsize = analysis->DofsPerNode(&this->doftype,iomodel->domaintype,in_approximation);
     63        delete analysis;
     64
    6965        if(this->gsize>0){
    70                 this->f_set    = xNew<bool>(this->gsize);
    71                 this->s_set    = xNew<bool>(this->gsize);
    72                 this->svalues  = xNew<IssmDouble>(this->gsize);
    73                 this->gdoflist = xNew<int>(this->gsize);
     66                this->f_set          = xNew<bool>(this->gsize);
     67                this->s_set          = xNew<bool>(this->gsize);
     68                this->svalues        = xNew<IssmDouble>(this->gsize);
     69                this->gdoflist       = xNew<int>(this->gsize);
     70                this->gdoflist_local = xNew<int>(this->gsize);
     71                this->fsize          = -1;
     72                this->ssize          = -1;
     73                this->fdoflist       = NULL;
     74                this->sdoflist       = NULL;
     75                this->fdoflist_local = NULL;
     76                this->sdoflist_local = NULL;
     77        }
     78        else{
     79                this->f_set          = NULL;
     80                this->s_set          = NULL;
     81                this->svalues        = NULL;
     82                this->gdoflist       = NULL;
     83                this->gdoflist_local = NULL;
     84                this->fsize          = -1;
     85                this->ssize          = -1;
     86                this->fdoflist       = NULL;
     87                this->sdoflist       = NULL;
     88                this->fdoflist_local = NULL;
     89                this->sdoflist_local = NULL;
    7490        }
    7591
     
    8197                this->gdoflist[i] = -1;
    8298        }
    83         delete analysis;
    84 
    85         if(analysis_enum==StressbalanceAnalysisEnum)
    86          this->approximation=in_approximation;
    87         else
    88          this->approximation=0;
    8999
    90100        /*Stop here if AMR*/
     
    136146                                analysis_enum==HydrologyDCEfficientAnalysisEnum ||
    137147                                analysis_enum==LevelsetAnalysisEnum
    138                                 ){
     148          ){
    139149                if(iomodel->domaintype!=Domain2DhorizontalEnum & iomodel->domaintype!=Domain3DsurfaceEnum){
    140150                        /*On a 3d mesh, we may have collapsed elements, hence dead nodes. Freeze them out: */
     
    147157        if(
    148158                                analysis_enum==FreeSurfaceTopAnalysisEnum
    149                                 ){
     159          ){
    150160                if(iomodel->domaintype!=Domain2DhorizontalEnum){
    151161                        /*On a 3d mesh, we may have collapsed elements, hence dead nodes. Freeze them out: */
     
    161171Node::~Node(){/*{{{*/
    162172
    163         if(this->f_set) xDelete<bool>(f_set);
    164         if(this->s_set) xDelete<bool>(s_set);
    165         if(this->svalues) xDelete<IssmDouble>(svalues);
    166         if(this->doftype) xDelete<int>(doftype);
    167         if(this->gdoflist) xDelete<int>(gdoflist);
    168         if(this->fdoflist) xDelete<int>(fdoflist);
    169         if(this->sdoflist) xDelete<int>(sdoflist);
     173        if(this->f_set)          xDelete<bool>(f_set);
     174        if(this->s_set)          xDelete<bool>(s_set);
     175        if(this->svalues)        xDelete<IssmDouble>(svalues);
     176        if(this->doftype)        xDelete<int>(doftype);
     177        if(this->gdoflist)       xDelete<int>(gdoflist);
     178        if(this->fdoflist)       xDelete<int>(fdoflist);
     179        if(this->sdoflist)       xDelete<int>(sdoflist);
     180        if(this->gdoflist_local) xDelete<int>(gdoflist_local);
     181        if(this->fdoflist_local) xDelete<int>(fdoflist_local);
     182        if(this->sdoflist_local) xDelete<int>(sdoflist_local);
    170183        return;
    171184}
     
    204217                if(this->doftype) output->doftype=xNew<int>(output->gsize);
    205218                output->gdoflist=xNew<int>(output->gsize);
    206         }
    207         if(output->fsize>0)output->fdoflist=xNew<int>(output->fsize);
    208         if(output->ssize>0)output->sdoflist=xNew<int>(output->ssize);
     219                output->gdoflist_local=xNew<int>(output->gsize);
     220        }
     221        if(output->fsize>0){
     222                output->fdoflist=xNew<int>(output->fsize);
     223                output->fdoflist_local=xNew<int>(output->fsize);
     224        }
     225        if(output->ssize>0){
     226                output->sdoflist=xNew<int>(output->ssize);
     227                output->sdoflist_local=xNew<int>(output->ssize);
     228        }
    209229
    210230        if(output->gsize>0){
     
    214234                if(output->doftype)memcpy(output->doftype,this->doftype,output->gsize*sizeof(int));
    215235                memcpy(output->gdoflist,this->gdoflist,output->gsize*sizeof(int));
    216         }
    217         if(output->fsize>0)memcpy(output->fdoflist,this->fdoflist,output->fsize*sizeof(int));
    218         if(output->ssize>0)memcpy(output->sdoflist,this->sdoflist,output->ssize*sizeof(int));
     236                memcpy(output->gdoflist_local,this->gdoflist_local,output->gsize*sizeof(int));
     237        }
     238        if(output->fsize>0){
     239                memcpy(output->fdoflist,this->fdoflist,output->fsize*sizeof(int));
     240                memcpy(output->fdoflist_local,this->fdoflist_local,output->fsize*sizeof(int));
     241        }
     242        if(output->ssize>0){
     243                memcpy(output->sdoflist,this->sdoflist,output->ssize*sizeof(int));
     244                memcpy(output->sdoflist_local,this->sdoflist_local,output->ssize*sizeof(int));
     245        }
    219246
    220247        return (Object*)output;
     
    245272        MARSHALLING_DYNAMIC(fdoflist,int,fsize);
    246273        MARSHALLING_DYNAMIC(sdoflist,int,ssize);
     274        MARSHALLING_DYNAMIC(gdoflist_local,int,gsize);
     275        MARSHALLING_DYNAMIC(fdoflist_local,int,fsize);
     276        MARSHALLING_DYNAMIC(sdoflist_local,int,ssize);
    247277} /*}}}*/
    248278
     
    284314
    285315        _printf_("   g_doflist (" << this->gsize << "): |");
    286         for(i=0;i<this->gsize;i++){
    287                 _printf_(" " << gdoflist[i] << " |");
    288         }
     316        for(i=0;i<this->gsize;i++) _printf_(" " << gdoflist[i] << " |");
    289317        _printf_("\n");
     318        _printf_("   g_doflist_local (" << this->gsize << "): |");
     319        for(i=0;i<this->gsize;i++) _printf_(" " << gdoflist_local[i] << " |");
     320        _printf_("\n");
    290321
    291322        _printf_("   f_doflist (" << this->fsize << "): |");
    292         for(i=0;i<this->fsize;i++){
    293                 _printf_(" " << fdoflist[i] << " |");
    294         }
     323        for(i=0;i<this->fsize;i++) _printf_(" " << fdoflist[i] << " |");
    295324        _printf_("\n");
     325        _printf_("   f_doflist_local (" << this->fsize << "): |");
     326        for(i=0;i<this->fsize;i++) _printf_(" " << fdoflist_local[i] << " |");
     327        _printf_("\n");
    296328
    297329        _printf_("   s_doflist (" << this->ssize << "): |");
    298         for(i=0;i<this->ssize;i++){
    299                 _printf_(" " << sdoflist[i] << " |");
    300         }
     330        for(i=0;i<this->ssize;i++) _printf_(" " << sdoflist[i] << " |");
     331        _printf_("\n");
     332        _printf_("   s_doflist_local (" << this->ssize << "): |");
     333        for(i=0;i<this->ssize;i++) _printf_(" " << sdoflist_local[i] << " |");
    301334        _printf_("\n");
    302335
     
    826859
    827860/* indexing routines:*/
    828 void Node::DistributeDofs(int* pdofcount,int setenum){/*{{{*/
    829 
    830         int i;
    831         int dofcount;
    832 
    833         dofcount=*pdofcount;
     861void Node::AllocateDofLists(int setenum){/*{{{*/
    834862
    835863        /*Initialize: */
    836864        int size=0;
     865
    837866        if(setenum==FsetEnum){
    838                 for(i=0;i<this->gsize;i++) if(f_set[i])size++;
     867                for(int i=0;i<this->gsize;i++) if(f_set[i])size++;
    839868                this->fsize=size;
    840869                xDelete<int>(this->fdoflist);
    841 
    842                 if(this->fsize)
    843                  this->fdoflist=xNew<int>(size);
    844                 else
    845                  this->fdoflist=NULL;
    846         }
    847         if(setenum==SsetEnum){
    848                 for(i=0;i<this->gsize;i++) if(s_set[i])size++;
     870                xDelete<int>(this->fdoflist_local);
     871
     872                if(this->fsize){
     873                        this->fdoflist=xNew<int>(size);
     874                        this->fdoflist_local=xNew<int>(size);
     875                }
     876                else{
     877                        this->fdoflist=NULL;
     878                        this->fdoflist_local=NULL;
     879                }
     880        }
     881        else if(setenum==SsetEnum){
     882                for(int i=0;i<this->gsize;i++) if(s_set[i])size++;
    849883                this->ssize=size;
    850884                xDelete<int>(this->sdoflist);
    851 
    852                 if(this->ssize)
    853                  this->sdoflist=xNew<int>(size);
    854                 else
    855                  this->sdoflist=NULL;
    856         }
    857 
    858         /*For clone nodfs, don't distribute dofs, we will get them from another cpu in UpdateCloneDofs!*/
    859         if(clone) return;
     885                xDelete<int>(this->sdoflist_local);
     886
     887                if(this->ssize){
     888                        this->sdoflist=xNew<int>(size);
     889                        this->sdoflist_local=xNew<int>(size);
     890                }
     891                else{
     892                        this->sdoflist=NULL;
     893                        this->sdoflist_local=NULL;
     894                }
     895        }
     896        /*TODO, we should never be here for GSET! add an assert and track down whether this happens*/
     897}
     898/*}}}*/
     899void Node::DistributeLocalDofs(int* pdofcount,int setenum){/*{{{*/
     900
     901        /*Get current count*/
     902        int dofcount=*pdofcount;
    860903
    861904        /*This node should distribute dofs for setenum set (eg, f_set or s_set), go ahead: */
    862905        if(setenum==GsetEnum){
    863                 for(i=0;i<this->gsize;i++){
    864                         gdoflist[i]=dofcount+i;
    865                 }
     906                _assert_(this->gsize==0 || this->gdoflist_local);
     907                for(int i=0;i<this->gsize;i++) gdoflist_local[i]=dofcount+i;
    866908                dofcount+=this->gsize;
    867909        }
    868910        else if(setenum==FsetEnum){
    869                 for(i=0;i<this->fsize;i++){
    870                         fdoflist[i]=dofcount+i;
    871                 }
     911                _assert_(this->fsize==0 || this->fdoflist_local);
     912                for(int i=0;i<this->fsize;i++) fdoflist_local[i]=dofcount+i;
    872913                dofcount+=this->fsize;
    873914        }
    874915        else if(setenum==SsetEnum){
    875                 for(i=0;i<this->ssize;i++){
    876                         sdoflist[i]=dofcount+i;
    877                 }
     916                _assert_(this->ssize==0 || this->sdoflist_local);
     917                for(int i=0;i<this->ssize;i++) sdoflist_local[i]=dofcount+i;
    878918                dofcount+=this->ssize;
    879919        }
     
    884924}
    885925/*}}}*/
    886 void Node::OffsetDofs(int dofcount,int setenum){/*{{{*/
    887 
    888         int i;
    889 
    890         if(clone){
    891                 /*This node is a clone, don't off_set the dofs!: */
    892                 return;
    893         }
     926void Node::DistributeGlobalDofsMasters(int dofcount,int setenum){/*{{{*/
     927
     928        /*This node is a clone, don't offset the dofs!: */
     929        if(clone) return;
    894930
    895931        /*This node should off_set the dofs, go ahead: */
    896932        if(setenum==GsetEnum){
    897                 for(i=0;i<this->gsize;i++) gdoflist[i]+=dofcount;
     933                _assert_(this->gsize==0 || this->gdoflist);
     934                for(int i=0;i<this->gsize;i++) this->gdoflist[i]=this->gdoflist_local[i]+dofcount;
    898935        }
    899936        else if(setenum==FsetEnum){
    900                 for(i=0;i<this->fsize;i++) fdoflist[i]+=dofcount;
     937                _assert_(this->fsize==0 || this->fdoflist);
     938                for(int i=0;i<this->fsize;i++) this->fdoflist[i]=this->fdoflist_local[i]+dofcount;
    901939        }
    902940        else if(setenum==SsetEnum){
    903                 for(i=0;i<this->ssize;i++) sdoflist[i]+=dofcount;
     941                _assert_(this->ssize==0 || this->sdoflist);
     942                for(int i=0;i<this->ssize;i++) this->sdoflist[i]=this->sdoflist_local[i]+dofcount;
    904943        }
    905944        else _error_("set of enum type " << EnumToStringx(setenum) << " not supported yet!");
    906945}
    907946/*}}}*/
    908 void Node::ShowTrueDofs(int* truedofs,int setenum){/*{{{*/
     947void Node::ShowMasterDofs(int* truedofs,int setenum){/*{{{*/
    909948
    910949        _assert_(!this->clone);
  • issm/trunk-jpl/src/c/classes/Node.h

    r23612 r23627  
    2424
    2525        private:
    26                 int approximation; //For ice flow models, we need to know what ice flow approximation is employed on this node
     26                int  approximation; //For ice flow models, we need to know what ice flow approximation is employed on this node
     27                bool clone;  //this node is replicated from another one
     28
     29                int  id;    // unique arbitrary id.
     30                int  sid;   // "serial" id (rank of this node if the dataset was serial on 1 cpu)
     31                int  lid;   // "local"  id (rank of this node in current partition)
     32                int  pid;   // parallel id (specific to this partition)
    2733
    2834        public:
    29 
    30                 int id;    // unique arbitrary id.
    31                 int sid;   // "serial" id (rank of this node if the dataset was serial on 1 cpu)
    32                 int lid;   // "local"  id (rank of this node in current partition)
    33                 int pid;   // parallel id (specific to this partition)
    34 
    35                 bool clone;  //this node is replicated from another one
    3635                int  analysis_enum;
    3736                IssmDouble   coord_system[3][3];
     
    5655
    5756                /*list of degrees of freedom: */
    58                 int *gdoflist;   //dof list in g_set
    59                 int *fdoflist;   //dof list in f_set
    60                 int *sdoflist;   //dof list in s_set
     57                int *gdoflist;
     58                int *fdoflist;
     59                int *sdoflist;
     60                int *gdoflist_local;
     61                int *fdoflist_local;
     62                int *sdoflist_local;
    6163
    6264                /*Node constructors, destructors*/
     
    7476
    7577                /*Node numerical routines*/
     78                void  AllocateDofLists(int setenum);
    7679                void  Activate(void);
    7780                void  ApplyConstraint(int dof,IssmDouble value);
    7881                void  CreateNodalConstraints(Vector<IssmDouble>* ys);
    7982                void  Deactivate(void);
    80                 void  DistributeDofs(int* pdofcount,int setenum);
     83                void  DistributeLocalDofs(int* pdofcount,int setenum);
    8184                void  DofInFSet(int dof);
    8285                void  DofInSSet(int dof);
     
    9497                int   IsGrounded();
    9598                int   Lid(void);
    96                 void  OffsetDofs(int dofcount,int setenum);
     99                void  DistributeGlobalDofsMasters(int dofcount,int setenum);
    97100                void  ReindexingDone(void);
    98101                void  RelaxConstraint(int dof);
    99102                bool  RequiresDofReindexing(void);
    100103                void  SetCurrentConfiguration(DataSet* nodes,Vertices* vertices);
    101                 void  ShowTrueDofs(int* truerows,int setenum);
     104                void  ShowMasterDofs(int* truerows,int setenum);
    102105                int   Sid(void);
    103106                void  UpdateCloneDofs(int* alltruerows,int setenum);
  • issm/trunk-jpl/src/c/classes/Nodes.cpp

    r23600 r23627  
    138138        int num_procs = IssmComm::GetSize();
    139139
    140         /*Go through objects, and distribute dofs locally, from 0 to numberofdofsperobject*/
     140        /*First, allocate dof lists based on fset and sset*/
     141        for(int i=0;i<this->Size();i++){
     142                Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i));
     143                node->AllocateDofLists(setenum);
     144        }
     145
     146        /*Now, Build local dofs for masters first*/
    141147        int  dofcount=0;
    142148        for(int i=0;i<this->Size();i++){
    143149                Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i));
    144                 node->DistributeDofs(&dofcount,setenum);
     150                if(!node->IsClone()) node->DistributeLocalDofs(&dofcount,setenum);
     151        }
     152        /*Build local dofs for clones, they always will be at the end*/
     153        int dofcount_local = dofcount;
     154        for(int i=0;i<this->Size();i++){
     155                Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i));
     156                if(node->IsClone()) node->DistributeLocalDofs(&dofcount_local,setenum);
    145157        }
    146158
    147159        /* Now every object has distributed dofs, but locally, and with a dof count starting from
    148160         * 0. This means the dofs between all the cpus are not unique. We now offset the dofs of each
    149          * cpus by the total last dofs of the previus cpu, starting from 0.
     161         * cpus by the total last (master) dofs of the previus cpu, starting from 0.
    150162         * First: get number of dofs for each cpu*/
    151163        int* alldofcount=xNew<int>(num_procs);
     
    161173                /*Check that this node corresponds to our analysis currently being carried out: */
    162174                Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(i));
    163                 node->OffsetDofs(offset,setenum);
     175                node->DistributeGlobalDofsMasters(offset,setenum);
    164176        }
    165177
     
    174186                        for(int i=0;i<numids;i++){
    175187                                Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(this->common_send_ids[rank][i]));
    176                                 node->ShowTrueDofs(&truedofs[i*maxdofspernode+0],setenum);
     188                                node->ShowMasterDofs(&truedofs[i*maxdofspernode+0],setenum);
    177189                        }
    178190                        ISSM_MPI_Send(truedofs,numids*maxdofspernode,ISSM_MPI_INT,rank,0,IssmComm::GetComm());
Note: See TracChangeset for help on using the changeset viewer.