Changeset 22491


Ignore:
Timestamp:
03/01/18 15:40:58 (7 years ago)
Author:
seroussi
Message:

ADD: starting to exchange info with MITgcm, will clean up later

File:
1 edited

Legend:

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

    r22471 r22491  
    2121        /*parameters: */
    2222        IssmDouble finaltime,dt,yts;
    23         bool       isstressbalance,ismasstransport,issmb,isFS,isthermal,isgroundingline,isgia,isslr,iscoupler,ismovingfront,isdamageevolution,ishydrology;
     23        bool       isstressbalance,ismasstransport,issmb,isFS,isthermal,isgroundingline,isgia,isslr,iscoupler,ismovingfront,isdamageevolution,ishydrology,isoceancoupling;
    2424        bool       save_results,dakota_analysis;
    2525        int        timestepping;
     
    5757        femmodel->parameters->FindParam(&isgroundingline,TransientIsgroundinglineEnum);
    5858        femmodel->parameters->FindParam(&ismovingfront,TransientIsmovingfrontEnum);
     59        femmodel->parameters->FindParam(&isoceancoupling,TransientIsoceancouplingEnum);
    5960        femmodel->parameters->FindParam(&isdamageevolution,TransientIsdamageevolutionEnum);
    6061        femmodel->parameters->FindParam(&ishydrology,TransientIshydrologyEnum);
     
    7172        }
    7273        #endif
     74
     75        if(isoceancoupling){
     76                printf("startin ocean coupling \n");
     77                int my_rank;
     78                int oceannxsize,oceannysize;
     79                IssmDouble *oceangridx;
     80                IssmDouble *oceangridy;
     81                IssmDouble *oceanmelt;
     82                IssmDouble *icebase;
     83                IssmDouble coupling_time,oceantime;
     84                ISSM_MPI_Comm tomitgcmcomm;
     85                ISSM_MPI_Status status;
     86
     87                my_rank=IssmComm::GetRank();
     88                femmodel->parameters->FindParam(&coupling_time,TimesteppingCouplingTimeEnum);
     89                GenericParam<ISSM_MPI_Comm>* parcom = dynamic_cast<GenericParam<ISSM_MPI_Comm>*>(femmodel->parameters->FindParamObject(ToMITgcmCommEnum));
     90                if(!parcom)_error_("TransferForcing error message: could not find ToMITgcmCommEnum communicator");
     91                tomitgcmcomm=parcom->GetParameterValue();
     92                if(my_rank==0){
     93                        printf("starting ocean coupling rank 0 \n");
     94                        printf("coupling time: %g \n",coupling_time);
     95                        ISSM_MPI_Send(&coupling_time,1,ISSM_MPI_DOUBLE,0,10001000,tomitgcmcomm);
     96                        ISSM_MPI_Send(&time,1,ISSM_MPI_DOUBLE,0,10001001,tomitgcmcomm);
     97                        ISSM_MPI_Recv(&oceantime,1,ISSM_MPI_DOUBLE,0,10001002,tomitgcmcomm,&status);
     98                        ISSM_MPI_Recv(&oceannxsize,1,ISSM_MPI_INT,0,10001003,tomitgcmcomm,&status);
     99                        ISSM_MPI_Recv(&oceannysize,1,ISSM_MPI_INT,0,10001004,tomitgcmcomm,&status);
     100                        //oceangridx = xNew<IssmDouble>(oceannxsize*oceannysize);
     101                        //ISSM_MPI_Recv(&oceangridx,oceannxsize*oceannysize,ISSM_MPI_DOUBLE,0,10001005,tomitgcmcomm,&status);
     102                        //oceangridy = xNew<IssmDouble>(oceannxsize*oceannysize);
     103                        //ISSM_MPI_Recv(&oceangridy,oceannxsize*oceannysize,ISSM_MPI_DOUBLE,0,10001006,tomitgcmcomm,&status);
     104                        //icebase = xNew<IssmDouble>(oceannxsize*oceannysize);
     105                        //ISSM_MPI_Send(&icebase,oceannxsize*oceannysize,ISSM_MPI_DOUBLE,0,10001008,tomitgcmcomm);
     106                        //oceanmelt = xNew<IssmDouble>(oceannxsize*oceannysize);
     107                        //ISSM_MPI_Recv(&oceanmelt,oceannxsize*oceannysize,ISSM_MPI_DOUBLE,0,10001007,tomitgcmcomm,&status);
     108                        //xDelete<IssmDouble>(icebase);
     109                        //xDelete<IssmDouble>(oceangridx);
     110                        //xDelete<IssmDouble>(oceangridy);
     111                        //xDelete<IssmDouble>(oceanmelt);
     112                }
     113        }
    73114
    74115        while(time < finaltime - (yts*DBL_EPSILON)){ //make sure we run up to finaltime.
Note: See TracChangeset for help on using the changeset viewer.