Changeset 27814


Ignore:
Timestamp:
06/30/23 14:01:03 (21 months ago)
Author:
seroussi
Message:

CHG: exchanging ice mass instead of ice base elevation for ocean model

File:
1 edited

Legend:

Unmodified
Added
Removed
  • TabularUnified issm/trunk-jpl/src/c/modules/OceanExchangeDatax/OceanExchangeDatax.cpp

    r27812 r27814  
    2828        IssmDouble *oceangridx;
    2929        IssmDouble *oceangridy;
    30         IssmDouble *icebase_oceangrid = NULL;
     30        IssmDouble *icethickness_oceangrid = NULL;
    3131        IssmDouble *icemask_oceangrid = NULL;
    3232        IssmDouble* x_ice             = NULL;
     
    3434        IssmDouble* lat_ice           = NULL;
    3535        IssmDouble* lon_ice           = NULL;
    36         IssmDouble* icebase           = NULL;
     36        IssmDouble* icethickness      = NULL;
    3737        IssmDouble* icemask           = NULL;
    3838        IssmDouble* melt_mesh         = NULL;
     
    8181        }
    8282
    83         /*Interpolate ice base and mask onto ocean grid*/
     83        /*Interpolate ice thickness and mask onto ocean grid*/
    8484        femmodel->GetMesh(femmodel->vertices,femmodel->elements,&x_ice,&y_ice,&index_ice);
    8585        BamgTriangulatex(&index_ocean,&nels_ocean,oceangridx,oceangridy,ngrids_ocean);
    8686        femmodel->vertices->XYList(&lon_ice,&lat_ice);
    8787        //femmodel->vertices->LatLonList(&lat_ice,&lon_ice);
    88         GetVectorFromInputsx(&icebase,femmodel,BaseEnum,VertexSIdEnum);
     88        GetVectorFromInputsx(&icethickness,femmodel,ThicknessEnum,VertexSIdEnum);
    8989        Options* options = new Options();
    9090        GenericOption<double> *odouble = new GenericOption<double>();
     
    9696        odouble->size[1]=1;
    9797        options->AddOption(odouble);
    98         InterpFromMeshToMesh2dx(&icebase_oceangrid,index_ice,lon_ice,lat_ice,ngrids_ice,nels_ice,
    99                                         icebase,ngrids_ice,1,oceangridx,oceangridy,ngrids_ocean,options);
     98        InterpFromMeshToMesh2dx(&icethickness_oceangrid,index_ice,lon_ice,lat_ice,ngrids_ice,nels_ice,
     99                                        icethickness,ngrids_ice,1,oceangridx,oceangridy,ngrids_ocean,options);
    100100        delete options;
    101         xDelete<IssmDouble>(icebase);
     101        xDelete<IssmDouble>(icethickness);
    102102
    103103        GetVectorFromInputsx(&icemask,femmodel,MaskIceLevelsetEnum,VertexSIdEnum);
     
    117117
    118118        /*Put +9999 for places where there is no ice!*/
    119         for(int i=0;i<ngrids_ocean;i++) if(icemask_oceangrid[i]>0.) icebase_oceangrid[i]=+9999.;
     119        femmodel->parameters->FindParam(&rho_ice,MaterialsRhoIceEnum);
     120        for(int i=0;i<ngrids_ocean;i++) icethickness_oceangrid[i]=icethickness_oceangrid[i]*rho_ice; //ocean needs ice mass in kg/m^2
     121        for(int i=0;i<ngrids_ocean;i++) if(icemask_oceangrid[i]>0.) icethickness_oceangrid[i]=+9999.;
    120122        xDelete<IssmDouble>(icemask_oceangrid);
    121123
    122         if(init_stage==true){ //just send icebase
     124        if(init_stage==true){ //just send icethickness
    123125                if(my_rank==0){
    124                         ISSM_MPI_Send(icebase_oceangrid,ngrids_ocean,ISSM_MPI_DOUBLE,0,10001008,tomitgcmcomm);
     126                        ISSM_MPI_Send(icethickness_oceangrid,ngrids_ocean,ISSM_MPI_DOUBLE,0,10001008,tomitgcmcomm);
    125127                }
    126128        }
    127129        else{ //send and receive exchanged data
    128                 femmodel->parameters->FindParam(&rho_ice,MaterialsRhoIceEnum);
    129130                femmodel->parameters->FindParam(&yts,ConstantsYtsEnum);
    130131                if(my_rank==0){
     
    134135                        oceanmelt = xNew<IssmDouble>(ngrids_ocean);
    135136                        ISSM_MPI_Recv(oceanmelt,ngrids_ocean,ISSM_MPI_DOUBLE,0,10001007,tomitgcmcomm,&status);
    136                         ISSM_MPI_Send(icebase_oceangrid,ngrids_ocean,ISSM_MPI_DOUBLE,0,10001008,tomitgcmcomm);
     137                        ISSM_MPI_Send(icethickness_oceangrid,ngrids_ocean,ISSM_MPI_DOUBLE,0,10001008,tomitgcmcomm);
    137138                }
    138139                ISSM_MPI_Bcast(&oceantime,1,ISSM_MPI_DOUBLE,0,IssmComm::GetComm());
     
    156157        xDelete<IssmDouble>(x_ice);
    157158        xDelete<IssmDouble>(y_ice);
    158         xDelete<IssmDouble>(icebase_oceangrid);
     159        xDelete<IssmDouble>(icethickness_oceangrid);
    159160        xDelete<IssmDouble>(oceangridx);
    160161        xDelete<IssmDouble>(oceangridy);
Note: See TracChangeset for help on using the changeset viewer.