Changeset 22709


Ignore:
Timestamp:
04/25/18 10:06:20 (7 years ago)
Author:
seroussi
Message:

CHG: changed the order of the coupling

File:
1 edited

Legend:

Unmodified
Added
Removed
  • issm/trunk-jpl/src/c/cores/transient_core.cpp

    r22703 r22709  
    8282                IssmDouble *oceangridx;
    8383                IssmDouble *oceangridy;
    84                 IssmDouble *oceanmelt;
    8584                IssmDouble *icebase;
    8685                IssmDouble coupling_time,oceantime;
     
    9594                if(my_rank==0){
    9695                        ISSM_MPI_Send(&coupling_time,1,ISSM_MPI_DOUBLE,0,10001000,tomitgcmcomm);
    97                         ISSM_MPI_Send(&time,1,ISSM_MPI_DOUBLE,0,10001001,tomitgcmcomm);
    98                         ISSM_MPI_Recv(&oceantime,1,ISSM_MPI_DOUBLE,0,10001002,tomitgcmcomm,&status);
    9996                        ISSM_MPI_Recv(&oceangridnxsize,1,ISSM_MPI_INT,0,10001003,tomitgcmcomm,&status);
    10097                        ISSM_MPI_Recv(&oceangridnysize,1,ISSM_MPI_INT,0,10001004,tomitgcmcomm,&status);
     
    105102                        oceangridy = xNew<IssmDouble>(oceangridnxsize*oceangridnysize);
    106103                        ISSM_MPI_Recv(oceangridy,oceangridnxsize*oceangridnysize,ISSM_MPI_DOUBLE,0,10001006,tomitgcmcomm,&status);
     104
     105                        ISSM_MPI_Send(&time,1,ISSM_MPI_DOUBLE,0,10001001,tomitgcmcomm);
     106                        ISSM_MPI_Recv(&oceantime,1,ISSM_MPI_DOUBLE,0,10001002,tomitgcmcomm,&status);
    107107                        icebase = xNew<IssmDouble>(oceangridnxsize*oceangridnysize);
    108108                        for(int i=0;i<oceangridnxsize;i++){
     
    112112                        }
    113113                        ISSM_MPI_Send(icebase,oceangridnxsize*oceangridnysize,ISSM_MPI_DOUBLE,0,10001008,tomitgcmcomm);
    114                         oceanmelt = xNew<IssmDouble>(oceangridnxsize*oceangridnysize);
    115                         ISSM_MPI_Recv(oceanmelt,oceangridnxsize*oceangridnysize,ISSM_MPI_DOUBLE,0,10001007,tomitgcmcomm,&status);
    116114                        xDelete<IssmDouble>(icebase);
    117115                        xDelete<IssmDouble>(oceangridx);
    118116                        xDelete<IssmDouble>(oceangridy);
    119                         xDelete<IssmDouble>(oceanmelt);
    120117                }
    121118        }
     
    166163
    167164                if(isoceancoupling){ {{{
    168                         if(VerboseSolution()) _printf0_("   ocean coupling: sending ice base\n");
     165                        if(VerboseSolution()) _printf0_("   ocean coupling: exchanging information\n");
    169166                        int my_rank;
    170167                        int oceangridnxsize,oceangridnysize;
     168                        IssmDouble *oceanmelt;
    171169                        IssmDouble *icebase;
    172170                        IssmDouble oceantime;
     
    184182                                femmodel->parameters->FindParam(&oceangridnxsize,OceanGridNxEnum);
    185183                                femmodel->parameters->FindParam(&oceangridnysize,OceanGridNyEnum);
     184                                oceanmelt = xNew<IssmDouble>(oceangridnxsize*oceangridnysize);
     185                                ISSM_MPI_Recv(oceanmelt,oceangridnxsize*oceangridnysize,ISSM_MPI_DOUBLE,0,10001007,tomitgcmcomm,&status);
    186186                                for(int i=0;i<oceangridnxsize;i++){
    187187                                        for(int j=0;j<oceangridnysize;j++){
     
    190190                                }
    191191                                ISSM_MPI_Send(icebase,oceangridnxsize*oceangridnysize,ISSM_MPI_DOUBLE,0,10001008,tomitgcmcomm);
     192                                xDelete<IssmDouble>(oceanmelt);
    192193                                xDelete<IssmDouble>(icebase);
    193194                        }
     
    272273                }
    273274
    274                 if(isoceancoupling){ {{{
    275                         if(VerboseSolution()) _printf0_("   ocean coupling: receiving melt\n");
    276                         int my_rank;
    277                         int oceangridnxsize,oceangridnysize;
    278                         IssmDouble *oceanmelt;
    279                         ISSM_MPI_Comm tomitgcmcomm;
    280                         ISSM_MPI_Status status;
    281 
    282                         my_rank=IssmComm::GetRank();
    283                         GenericParam<ISSM_MPI_Comm>* parcom = dynamic_cast<GenericParam<ISSM_MPI_Comm>*>(femmodel->parameters->FindParamObject(ToMITgcmCommEnum));
    284                         if(!parcom)_error_("TransferForcing error message: could not find ToMITgcmCommEnum communicator");
    285                         tomitgcmcomm=parcom->GetParameterValue();
    286                         if(my_rank==0){
    287                                 femmodel->parameters->FindParam(&oceangridnxsize,OceanGridNxEnum);
    288                                 femmodel->parameters->FindParam(&oceangridnysize,OceanGridNyEnum);
    289                                 oceanmelt = xNew<IssmDouble>(oceangridnxsize*oceangridnysize);
    290                                 ISSM_MPI_Recv(oceanmelt,oceangridnxsize*oceangridnysize,ISSM_MPI_DOUBLE,0,10001007,tomitgcmcomm,&status);
    291                                 xDelete<IssmDouble>(oceanmelt);
    292                         }
    293                 }
    294                 }}}
    295 
    296275                if(recording_frequency && (step%recording_frequency==0)){
    297276                        if(VerboseSolution()) _printf0_("   checkpointing model \n");
Note: See TracChangeset for help on using the changeset viewer.