Changeset 23636


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

CHG: moving function to Femmodel so that it can be recycled

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

Legend:

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

    r23598 r23636  
    13061306
    13071307}/*}}}*/
     1308void FemModel::GetLocalVectorWithClonesGset(IssmDouble** plocal_ug,Vector<IssmDouble> *ug){/*{{{*/
     1309
     1310        /*recover my_rank:*/
     1311        ISSM_MPI_Status status;
     1312        int my_rank   = IssmComm::GetRank();
     1313        int num_procs = IssmComm::GetSize();
     1314
     1315        /*retrieve node info*/
     1316        int glocalsize         = this->nodes->NumberOfDofsLocalAll(GsetEnum);
     1317        int glocalsize_masters = this->nodes->NumberOfDofsLocal(GsetEnum);
     1318        int maxdofspernode     = this->nodes->MaxNumDofs(GsetEnum);
     1319
     1320        /*Get local vector of ug*/
     1321        int        *indices_ug_masters = NULL;
     1322        IssmDouble *local_ug_masters   = NULL;
     1323        ug->GetLocalVector(&local_ug_masters,&indices_ug_masters);
     1324        _assert_(glocalsize_masters==indices_ug_masters[glocalsize_masters-1] - indices_ug_masters[0]+1);
     1325        xDelete<int>(indices_ug_masters);
     1326
     1327        /*Now, extend vectors to account for clones (make vectors longer, for clones at the end)*/
     1328        IssmDouble *local_ug  = xNew<IssmDouble>(glocalsize);
     1329        xMemCpy<IssmDouble>(local_ug,local_ug_masters,glocalsize_masters);
     1330        xDelete<IssmDouble>(local_ug_masters);
     1331
     1332        /*Now send and receive ug for nodes on partition edge*/
     1333        #ifdef _HAVE_AD_
     1334        IssmDouble* buffer = xNew<IssmDouble>(this->nodes->Size()*maxdofspernode,"t"); //only one alloc, "t" is required by adolc
     1335        #else
     1336        IssmDouble* buffer = xNew<IssmDouble>(this->nodes->Size()*maxdofspernode);
     1337        #endif
     1338        for(int rank=0;rank<num_procs;rank++){
     1339                if(this->nodes->common_send[rank]){
     1340                        int  numids = this->nodes->common_send[rank];
     1341                        for(int i=0;i<numids;i++){
     1342                                int   master_lid = this->nodes->common_send_ids[rank][i];
     1343                                Node* node=xDynamicCast<Node*>(this->nodes->GetObjectByOffset(master_lid));
     1344                                _assert_(!node->IsClone());
     1345                                for(int j=0;j<node->gsize;j++) buffer[i*maxdofspernode+j]=local_ug[node->gdoflist_local[j]];
     1346                        }
     1347                        ISSM_MPI_Send(buffer,numids*maxdofspernode,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm());
     1348                }
     1349        }
     1350        for(int rank=0;rank<num_procs;rank++){
     1351                if(this->nodes->common_recv[rank]){
     1352                        int  numids = this->nodes->common_recv[rank];
     1353                        ISSM_MPI_Recv(buffer,numids*maxdofspernode,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
     1354                        for(int i=0;i<numids;i++){
     1355                                int   master_lid = this->nodes->common_recv_ids[rank][i];
     1356                                Node* node=xDynamicCast<Node*>(this->nodes->GetObjectByOffset(master_lid));
     1357                                for(int j=0;j<node->gsize;j++) local_ug[node->gdoflist_local[j]] = buffer[i*maxdofspernode+j];
     1358                        }
     1359                }
     1360        }
     1361        xDelete<IssmDouble>(buffer);
     1362
     1363        /*Assign output pointer*/
     1364        *plocal_ug = local_ug;
     1365}/*}}}*/
    13081366void FemModel::GroundedAreax(IssmDouble* pV, bool scaled){/*{{{*/
    13091367
  • issm/trunk-jpl/src/c/classes/FemModel.h

    r23585 r23636  
    9696                void FloatingAreax(IssmDouble* pV, bool scaled);
    9797                void GetInputLocalMinMaxOnNodesx(IssmDouble** pmin,IssmDouble** pmax,IssmDouble* ug);
     98                void GetLocalVectorWithClonesGset(IssmDouble** plocal_ug,Vector<IssmDouble> *ug);
    9899                void GroundedAreax(IssmDouble* pV, bool scaled);
    99100                void IceMassx(IssmDouble* pV, bool scaled);
  • issm/trunk-jpl/src/c/modules/InputUpdateFromSolutionx/InputUpdateFromSolutionx.cpp

    r23632 r23636  
    1717        Analysis* analysis = EnumToAnalysis(analysisenum);
    1818
    19         /*recover my_rank:*/
    20         ISSM_MPI_Status status;
    21         int my_rank   = IssmComm::GetRank();
    22         int num_procs = IssmComm::GetSize();
    23 
    24         /*retrieve node info*/
    25         int gsize              = femmodel->nodes->NumberOfDofs(GsetEnum);
    26         int glocalsize_masters = femmodel->nodes->NumberOfDofsLocal(GsetEnum);
    27         int glocalsize         = femmodel->nodes->NumberOfDofsLocalAll(GsetEnum);
    28         int maxdofspernode     = femmodel->nodes->MaxNumDofs(GsetEnum);
    29 
    30         /*Get local vector of ug*/
    31         int        *indices_ug_masters = NULL;
    32         IssmDouble *local_ug_masters   = NULL;
    33         solution->GetLocalVector(&local_ug_masters,&indices_ug_masters);
    34         _assert_(glocalsize_masters==indices_ug_masters[glocalsize_masters-1] - indices_ug_masters[0]+1);
    35         xDelete<int>(indices_ug_masters);
    36 
    37         /*Now, extend vectors to account for clones (make vectors longer, for clones at the end)*/
    38         IssmDouble *local_ug  = xNew<IssmDouble>(glocalsize);
    39         xMemCpy<IssmDouble>(local_ug,local_ug_masters,glocalsize_masters);
    40         xDelete<IssmDouble>(local_ug_masters);
    41 
    42         /*Now send and receive ug for nodes on partition edge*/
    43         #ifdef _HAVE_AD_
    44         IssmDouble* buffer = xNew<IssmDouble>(femmodel->nodes->Size()*maxdofspernode,"t"); //only one alloc, "t" is required by adolc
    45         #else
    46         IssmDouble* buffer = xNew<IssmDouble>(femmodel->nodes->Size()*maxdofspernode);
    47         #endif
    48         for(int rank=0;rank<num_procs;rank++){
    49                 if(femmodel->nodes->common_send[rank]){
    50                         int  numids = femmodel->nodes->common_send[rank];
    51                         for(int i=0;i<numids;i++){
    52                                 int   master_lid = femmodel->nodes->common_send_ids[rank][i];
    53                                 Node* node=xDynamicCast<Node*>(femmodel->nodes->GetObjectByOffset(master_lid));
    54                                 _assert_(!node->IsClone());
    55                                 for(int j=0;j<node->gsize;j++) buffer[i*maxdofspernode+j]=local_ug[node->gdoflist_local[j]];
    56                         }
    57                         ISSM_MPI_Send(buffer,numids*maxdofspernode,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm());
    58                 }
    59         }
    60         for(int rank=0;rank<num_procs;rank++){
    61                 if(femmodel->nodes->common_recv[rank]){
    62                         int  numids = femmodel->nodes->common_recv[rank];
    63                         ISSM_MPI_Recv(buffer,numids*maxdofspernode,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
    64                         for(int i=0;i<numids;i++){
    65                                 int   master_lid = femmodel->nodes->common_recv_ids[rank][i];
    66                                 Node* node=xDynamicCast<Node*>(femmodel->nodes->GetObjectByOffset(master_lid));
    67                                 for(int j=0;j<node->gsize;j++) local_ug[node->gdoflist_local[j]] = buffer[i*maxdofspernode+j];
    68                         }
    69                 }
    70         }
    71         xDelete<IssmDouble>(buffer);
     19        /*Get local vector with both masters and slaves:*/
     20        IssmDouble *local_ug = NULL;
     21        femmodel->GetLocalVectorWithClonesGset(&local_ug,solution);
    7222
    7323        /*Now update inputs (analysis specific)*/
  • issm/trunk-jpl/src/c/modules/SetActiveNodesLSMx/SetActiveNodesLSMx.cpp

    r23629 r23636  
    109109        vec_mask_ice->Assemble();
    110110
    111 
    112         /*FIXME: What follows should be copied and pasted into InputUpdateFromSolution Nodes*/
    113 
    114         /*recover my_rank:*/
    115         ISSM_MPI_Status status;
    116         int my_rank   = IssmComm::GetRank();
    117         int num_procs = IssmComm::GetSize();
    118 
    119         /*retrieve node info*/
    120         int glocalsize         = femmodel->nodes->NumberOfDofsLocalAll(GsetEnum);
    121         int maxdofspernode     = femmodel->nodes->MaxNumDofs(GsetEnum);
    122 
    123         /*Get local vector of ug*/
    124         int        *indices_ug_masters = NULL;
    125         IssmDouble *local_ug_masters   = NULL;
    126         vec_mask_ice->GetLocalVector(&local_ug_masters,&indices_ug_masters);
    127         _assert_(glocalsize_masters==indices_ug_masters[glocalsize_masters-1] - indices_ug_masters[0]+1);
    128         xDelete<int>(indices_ug_masters);
     111        /*Get local vector with masters and slaves*/
     112        IssmDouble *local_ug = NULL;
     113        femmodel->GetLocalVectorWithClonesGset(&local_ug,vec_mask_ice);
    129114        delete vec_mask_ice;
    130 
    131         /*Now, extend vectors to account for clones (make vectors longer, for clones at the end)*/
    132         IssmDouble *local_ug  = xNew<IssmDouble>(glocalsize);
    133         xMemCpy<IssmDouble>(local_ug,local_ug_masters,glocalsize_masters);
    134         xDelete<IssmDouble>(local_ug_masters);
    135 
    136         /*Now send and receive ug for nodes on partition edge*/
    137         IssmDouble* buffer = xNew<IssmDouble>(femmodel->nodes->Size()*maxdofspernode); //only one alloc
    138         for(int rank=0;rank<num_procs;rank++){
    139                 if(femmodel->nodes->common_send[rank]){
    140                         int  numids = femmodel->nodes->common_send[rank];
    141                         for(int i=0;i<numids;i++){
    142                                 int   master_lid = femmodel->nodes->common_send_ids[rank][i];
    143                                 Node* node=xDynamicCast<Node*>(femmodel->nodes->GetObjectByOffset(master_lid));
    144                                 _assert_(!node->IsClone());
    145                                 for(int j=0;j<node->gsize;j++) buffer[i*maxdofspernode+j]=local_ug[node->gdoflist_local[j]];
    146                         }
    147                         ISSM_MPI_Send(buffer,numids*maxdofspernode,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm());
    148                 }
    149         }
    150         for(int rank=0;rank<num_procs;rank++){
    151                 if(femmodel->nodes->common_recv[rank]){
    152                         int  numids = femmodel->nodes->common_recv[rank];
    153                         ISSM_MPI_Recv(buffer,numids*maxdofspernode,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
    154                         for(int i=0;i<numids;i++){
    155                                 int   master_lid = femmodel->nodes->common_recv_ids[rank][i];
    156                                 Node* node=xDynamicCast<Node*>(femmodel->nodes->GetObjectByOffset(master_lid));
    157                                 for(int j=0;j<node->gsize;j++) local_ug[node->gdoflist_local[j]] = buffer[i*maxdofspernode+j];
    158                         }
    159                 }
    160         }
    161         xDelete<IssmDouble>(buffer);
    162115
    163116        /*Now update inputs (analysis specific)*/
Note: See TracChangeset for help on using the changeset viewer.