Changeset 23521


Ignore:
Timestamp:
12/07/18 21:13:06 (6 years ago)
Author:
Mathieu Morlighem
Message:

BUG: fixing dakota, copy method of vertices was missing

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

Legend:

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

    r23520 r23521  
    3737        int num_proc=IssmComm::GetSize();
    3838
    39         if(common_recv) xDelete<int>(common_recv);
    40         if(common_send) xDelete<int>(common_send);
    41         if(common_recv_ids){
     39        if(this->common_recv){
     40                xDelete<int>(common_recv);
     41        }
     42        if(this->common_send) xDelete<int>(common_send);
     43        if(this->common_recv_ids){
    4244                for(int i=0;i<num_proc;i++) if(common_recv_ids[i]) xDelete<int>(common_recv_ids[i]);
    4345                xDelete<int*>(common_recv_ids);
    4446        }
    45         if(common_send_ids){
     47        if(this->common_send_ids){
    4648                for(int i=0;i<num_proc;i++) if(common_send_ids[i]) xDelete<int>(common_send_ids[i]);
    4749                xDelete<int*>(common_send_ids);
     
    5254
    5355/*Numerics management*/
     56Vertices* Vertices::Copy() {/*{{{*/
     57
     58        int num_proc = IssmComm::GetSize();
     59
     60        /*Copy dataset*/
     61        Vertices* output=new Vertices();
     62        output->sorted=this->sorted;
     63        output->numsorted=this->numsorted;
     64        output->presorted=this->presorted;
     65        for(vector<Object*>::iterator obj=this->objects.begin() ; obj < this->objects.end(); obj++ ){
     66                output->AddObject((*obj)->copy());
     67        }
     68
     69        /*Build id_offsets and sorted_ids*/
     70        int objsize = this->numsorted;
     71        if(this->sorted && objsize>0 && this->id_offsets){     
     72                /*Allocate new ids*/
     73                output->id_offsets=xNew<int>(objsize);
     74                xMemCpy<int>(output->id_offsets,this->id_offsets,objsize);
     75        }
     76        else output->id_offsets=NULL;
     77        if(this->sorted && objsize>0 && this->sorted_ids){
     78                /*Allocate new ids*/
     79                output->sorted_ids=xNew<int>(objsize);
     80                xMemCpy<int>(output->sorted_ids,this->sorted_ids,objsize);
     81        }
     82        else output->sorted_ids=NULL;
     83
     84
     85        if(this->common_recv){
     86                output->common_recv=xNew<int>(num_proc);
     87                for(int i=0;i<num_proc;i++) output->common_recv[i]=this->common_recv[i];
     88        }
     89        if(this->common_send){
     90                output->common_send=xNew<int>(num_proc);
     91                for(int i=0;i<num_proc;i++) output->common_send[i]=this->common_send[i];
     92        }
     93        if(this->common_recv_ids){
     94                output->common_recv_ids = xNew<int*>(num_proc);
     95                for(int i=0;i<num_proc;i++){
     96                        output->common_recv_ids[i]=xNew<int>(this->common_recv[i]);
     97                        for(int j=0;j<this->common_recv[i];j++) output->common_recv_ids[i][j]=this->common_recv_ids[i][j];
     98                }
     99        }
     100        if(this->common_send_ids){
     101                output->common_send_ids = xNew<int*>(num_proc);
     102                for(int i=0;i<num_proc;i++){
     103                        output->common_send_ids[i]=xNew<int>(this->common_send[i]);
     104                        for(int j=0;j<this->common_send[i];j++) output->common_send_ids[i][j]=this->common_send_ids[i][j];
     105                }
     106        }
     107
     108        return output;
     109}
     110/*}}}*/
    54111void Vertices::Marshall(char** pmarshalled_data,int* pmarshalled_data_size, int marshall_direction){ /*{{{*/
    55112
  • issm/trunk-jpl/src/c/classes/Vertices.h

    r23518 r23521  
    2525
    2626                /*Objects virtual functions*/
     27                Vertices* Copy();
    2728                void Marshall(char** pmarshalled_data,int* pmarshalled_data_size, int marshall_direction);
    2829
Note: See TracChangeset for help on using the changeset viewer.