Changeset 15861


Ignore:
Timestamp:
08/21/13 15:44:01 (12 years ago)
Author:
Mathieu Morlighem
Message:

NEW: added node->indexingupdate flag so that dofs are not renumbered except if absolutely necessary

Location:
issm/trunk-jpl/src/c
Files:
5 edited

Legend:

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

    r15771 r15861  
    3939
    4040        /*indexing:*/
     41        this->indexingupdate = true;
    4142        DistributeNumDofs(&this->indexing,analysis_type,in_approximation); //number of dofs per node
    4243        gsize=this->indexing.gsize;
     
    127128        _printf_("   analysis_type: " << EnumToStringx(analysis_type) << "\n");
    128129        _printf_("   approximation: " << EnumToStringx(approximation) << "\n");
     130        _printf_("   indexingupdate: " << indexingupdate << "\n");
    129131        indexing.Echo();
    130132        _printf_("   inputs: " << inputs << "\n");
     
    139141        _printf_("   sid: " << sid << "\n");
    140142        _printf_("   analysis_type: " << EnumToStringx(analysis_type) << "\n");
     143        _printf_("   approximation: " << EnumToStringx(approximation) << "\n");
     144        _printf_("   indexingupdate: " << indexingupdate << "\n");
    141145        indexing.DeepEcho();
    142146        _printf_("   inputs\n");
     
    159163int   Node::GetDof(int dofindex,int setenum){
    160164
     165        _assert_(!this->indexingupdate);
    161166        if(setenum==GsetEnum){
    162167                _assert_(dofindex>=0 && dofindex<indexing.gsize);
     
    179184        int count=0;
    180185        int count2=0;
     186
     187        _assert_(!this->indexingupdate);
    181188
    182189        if(approximation_enum==NoneApproximationEnum){
     
    242249        int count2=0;
    243250
     251        _assert_(!this->indexingupdate);
     252
    244253        if(approximation_enum==NoneApproximationEnum){
    245254                if(setenum==GsetEnum)for(i=0;i<this->indexing.gsize;i++) outdoflist[i]=i;
     
    365374}
    366375/*}}}*/
     376/*FUNCTION Node::RequiresDofReindexing{{{*/
     377bool Node::RequiresDofReindexing(void){
     378
     379        return this->indexingupdate;
     380
     381}
     382/*}}}*/
     383/*FUNCTION Node::ReindexingDone{{{*/
     384void Node::ReindexingDone(void){
     385
     386        this->indexingupdate = false;
     387
     388}
     389/*}}}*/
    367390/*FUNCTION Node::RelaxConstraint{{{*/
    368391void  Node::RelaxConstraint(int dof){
     
    376399void  Node::CreateVecSets(Vector<IssmDouble>* pv_g,Vector<IssmDouble>* pv_f,Vector<IssmDouble>* pv_s){
    377400
    378         IssmDouble gvalue=1.0; //all nodes are in the g set;
     401        IssmDouble gvalue=1.; //all nodes are in the g set;
    379402        IssmDouble value;
    380403
    381         int i;
    382 
    383         for(i=0;i<this->indexing.gsize;i++){
     404        for(int i=0;i<this->indexing.gsize;i++){
    384405
    385406                /*g set: */
     
    457478void  Node::Deactivate(void){
    458479
    459         indexing.Deactivate();
     480        if(IsActive()){
     481                this->indexingupdate = true;
     482                indexing.Deactivate();
     483        }
    460484
    461485}
     
    464488void  Node::Activate(void){
    465489
    466         indexing.Activate();
     490        if(!IsActive()){
     491                this->indexingupdate = true;
     492                indexing.Activate();
     493        }
    467494
    468495}
     
    491518
    492519        if(approximation_enum==NoneApproximationEnum){
    493                 if (setenum==GsetEnum) numdofs=this->indexing.gsize;
     520                if      (setenum==GsetEnum) numdofs=this->indexing.gsize;
    494521                else if (setenum==FsetEnum) numdofs=this->indexing.fsize;
    495522                else if (setenum==SsetEnum) numdofs=this->indexing.ssize;
     
    621648void   Node::UpdateSpcs(IssmDouble* ys){
    622649
    623         int     count=0;
    624         int     i;
    625 
    626         count=0;
    627         for(i=0;i<this->indexing.gsize;i++){
     650        int count = 0;
     651        for(int i=0;i<this->indexing.gsize;i++){
    628652                if(this->indexing.s_set[i]){
    629653                        this->indexing.svalues[i]=ys[this->indexing.sdoflist[count]];
     
    634658/*}}}*/
    635659/*FUNCTION Node::VecMerge {{{*/
    636 void   Node::VecMerge(Vector<IssmDouble>* ug, IssmDouble* vector_serial,int setenum){
    637 
    638         IssmDouble* values=NULL;
    639         int*    indices=NULL;
    640         int     count=0;
    641         int     i;
     660void Node::VecMerge(Vector<IssmDouble>* ug, IssmDouble* vector_serial,int setenum){
     661
     662        IssmDouble *values  = NULL;
     663        int        *indices = NULL;
     664        int         count   = 0;
     665        int         i;
    642666
    643667        if(setenum==FsetEnum){
  • issm/trunk-jpl/src/c/classes/Node.h

    r15771 r15861  
    3333                int sid;   //"serial" id (rank of this node if the dataset was serial on 1 cpu)
    3434
     35                bool         indexingupdate;
    3536                DofIndexing  indexing;
    3637                Inputs      *inputs;               //properties of this node
     
    8889                void  Deactivate(void);
    8990                void  UpdateSpcs(IssmDouble* ys);
     91                void  ReindexingDone(void);
     92                bool  RequiresDofReindexing(void);
    9093                int   IsFloating();
    9194                int   IsGrounded();
  • issm/trunk-jpl/src/c/classes/Nodes.cpp

    r15838 r15861  
    5151
    5252        /*Go through objects, and distribute dofs locally, from 0 to numberofdofsperobject*/
    53         for (i=0;i<this->Size();i++){
    54                 Node* node=dynamic_cast<Node*>(this->GetObjectByOffset(i));
    55 
    56                 /*Check that this node corresponds to our analysis currently being carried out: */
    57                 if (node->InAnalysis(analysis_type)){
     53        for(i=0;i<this->Size();i++){
     54                Node* node=dynamic_cast<Node*>(this->GetObjectByOffset(i));
     55
     56                /*Check that this node corresponds to our analysis currently being carried out: */
     57                if(node->InAnalysis(analysis_type)){
    5858                        node->DistributeDofs(&dofcount,setenum);
    5959                }
     
    7373                dofcount+=alldofcount[i];
    7474        }
    75         for (i=0;i<this->Size();i++){
     75        for(i=0;i<this->Size();i++){
    7676                /*Check that this node corresponds to our analysis currently being carried out: */
    7777                Node* node=dynamic_cast<Node*>(this->GetObjectByOffset(i));
     
    9191        }
    9292
    93         for (i=0;i<this->Size();i++){
     93        for(i=0;i<this->Size();i++){
    9494                Node* node=dynamic_cast<Node*>(this->GetObjectByOffset(i));
    9595                if (node->InAnalysis(analysis_type)){
     
    101101
    102102        /* Now every cpu knows the true dofs of everyone else that is not a clone*/
    103         for (i=0;i<this->Size();i++){
     103        for(i=0;i<this->Size();i++){
    104104                Node* node=dynamic_cast<Node*>(this->GetObjectByOffset(i));
    105105                if (node->InAnalysis(analysis_type)){
    106106                        node->UpdateCloneDofs(alltruedofs,maxdofspernode,setenum);
     107                }
     108        }
     109
     110        /*Update indexingupdateflag*/
     111        for(i=0;i<this->Size();i++){
     112                Node* node=dynamic_cast<Node*>(this->GetObjectByOffset(i));
     113                if (node->InAnalysis(analysis_type)){
     114                        node->ReindexingDone();
    107115                }
    108116        }
     
    330338}
    331339/*}}}*/
     340/*FUNCTION Nodes::RequiresDofReindexing{{{*/
     341bool Nodes::RequiresDofReindexing(int analysis_type){
     342
     343        int flag = 0;
     344        int allflag;
     345
     346        /*Now go through all nodes, and get how many dofs they own, unless they are clone nodes: */
     347        for(int i=0;i<this->Size();i++){
     348
     349                Node* node=dynamic_cast<Node*>(this->GetObjectByOffset(i));
     350
     351                /*Check that this node corresponds to our analysis currently being carried out: */
     352                if(node->InAnalysis(analysis_type)){
     353                        if(node->RequiresDofReindexing()){
     354                                flag = 1;
     355                                break;
     356                        }
     357                }
     358        }
     359
     360        /*Grab max of all cpus: */
     361        ISSM_MPI_Allreduce((void*)&flag,(void*)&allflag,1,ISSM_MPI_INT,ISSM_MPI_MAX,IssmComm::GetComm());
     362
     363        if(allflag){
     364                return true;
     365        }
     366        else{
     367                return false;
     368        }
     369}
     370/*}}}*/
  • issm/trunk-jpl/src/c/classes/Nodes.h

    r15634 r15861  
    2727                void  DistributeDofs(int analysis_type,int SETENUM);
    2828                void  FlagClones(int analysis_type);
     29                bool  RequiresDofReindexing(int analysis_type);
    2930                int   MaxNumDofs(int analysis_type,int setenum);
    3031                int   MaximumId(void);
  • issm/trunk-jpl/src/c/modules/NodesDofx/NodesDofx.cpp

    r14999 r15861  
    1010void NodesDofx(Nodes* nodes, Parameters* parameters,int configuration_type){
    1111
    12         int noerr=1;
    13         int found=0;
     12        /*Do we have any nodes for this analysis type? :*/
     13        if(!nodes->NumberOfNodes(configuration_type)) return;
    1414
    15         /*Do we have any nodes for this analysis type? :*/
    16         if(nodes->NumberOfNodes(configuration_type)){
     15        /*Do we really need to update dof indexings*/
     16        if(!nodes->RequiresDofReindexing(configuration_type)) return;
    1717
    18                 /*Ensure that only for each cpu, the partition border nodes only will be taken into account once
    19                  * across the cluster. To do so, we flag all the clone nodes: */
    20                 nodes->FlagClones(configuration_type);
     18        if(VerboseModule()) _printf0_("   Renumbering degrees of freedom\n");
    2119
    22                 /*Go through all nodes, and build degree of freedom lists. Each node gets a fixed number of dofs. When
    23                  *a  node has already been distributed dofs on one cpu, all other cpus with the same node cannot distribute it
    24                  *anymore. Use clone field to be sure of that: */
    25                 nodes->DistributeDofs(configuration_type,GsetEnum);
    26                 nodes->DistributeDofs(configuration_type,FsetEnum);
    27                 nodes->DistributeDofs(configuration_type,SsetEnum);
    28         }
     20        /*Ensure that only for each cpu, the partition border nodes only will be taken into account once
     21         * across the cluster. To do so, we flag all the clone nodes: */
     22        nodes->FlagClones(configuration_type);
     23
     24        /*Go through all nodes, and build degree of freedom lists. Each node gets a fixed number of dofs. When
     25         *a  node has already been distributed dofs on one cpu, all other cpus with the same node cannot distribute it
     26         *anymore. Use clone field to be sure of that: */
     27        nodes->DistributeDofs(configuration_type,GsetEnum);
     28        nodes->DistributeDofs(configuration_type,FsetEnum);
     29        nodes->DistributeDofs(configuration_type,SsetEnum);
    2930
    3031}
Note: See TracChangeset for help on using the changeset viewer.