Changeset 23080


Ignore:
Timestamp:
08/08/18 09:10:04 (7 years ago)
Author:
seroussi
Message:

CHG: added default value for ice base when no ice is present

File:
1 edited

Legend:

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

    r23079 r23080  
    9393                IssmDouble *oceangridx;
    9494                IssmDouble *oceangridy;
    95                 IssmDouble *icebase_oceangrid;
    96                 IssmDouble* x_ice = NULL;
    97                 IssmDouble* y_ice = NULL;
    98                 IssmDouble* lat_ice = NULL;
    99                 IssmDouble* lon_ice = NULL;
    100                 IssmDouble* z_ice = NULL;
    101                 IssmDouble* icebase= NULL;
    102                 int*        index_ice= NULL;
    103                 int*        index_ocean = NULL;
     95                IssmDouble *icebase_oceangrid = NULL;
     96                IssmDouble *icemask_oceangrid = NULL;
     97                IssmDouble* x_ice             = NULL;
     98                IssmDouble* y_ice             = NULL;
     99                IssmDouble* lat_ice           = NULL;
     100                IssmDouble* lon_ice           = NULL;
     101                IssmDouble* z_ice             = NULL;
     102                IssmDouble* icebase           = NULL;
     103                IssmDouble* icemask           = NULL;
     104                int*        index_ice         = NULL;
     105                int*        index_ocean       = NULL;
    104106                int         ngrids_ice=femmodel->vertices->NumberOfVertices();
    105107                int         nels_ice=femmodel->elements->NumberOfElements();
     
    138140                femmodel->parameters->SetParam(oceangridy,ngrids_ocean,OceanGridYEnum);
    139141
    140                 /*Interpolate ice base onto ocean grid*/
     142                /*Interpolate ice base and mask onto ocean grid*/
    141143                femmodel->GetMesh(femmodel->vertices,femmodel->elements,&x_ice,&y_ice,&z_ice,&index_ice);
    142144                BamgTriangulatex(&index_ocean,&nels_ocean,oceangridx,oceangridy,ngrids_ocean);
     
    157159                                                oceangridx,oceangridy,ngrids_ocean,options);
    158160                delete options;
    159 
     161                xDelete<IssmDouble>(icebase);
     162
     163                GetVectorFromInputsx(&icemask,femmodel,MaskIceLevelsetEnum,VertexSIdEnum);
     164                Options* options2 = new Options();
     165                GenericOption<double> *odouble2 = new GenericOption<double>();
     166                const char* name2 = "default";
     167                odouble2->name =xNew<char>(strlen(name2)+1);
     168                memcpy(odouble2->name,name2,(strlen(name2)+1)*sizeof(char));
     169                odouble2->value=+1.;
     170                odouble2->numel=1;
     171                odouble2->ndims=1;
     172                odouble2->size=NULL;
     173                options2->AddOption(odouble2);
     174                InterpFromMeshToMesh2dx(&icemask_oceangrid,index_ice,lon_ice,lat_ice,ngrids_ice,nels_ice,
     175                                        icemask,ngrids_ice,1,oceangridx,oceangridy,ngrids_ocean,options2);
     176                delete options2;
     177                xDelete<IssmDouble>(icemask);
     178
     179                /*Put +9999 for places where there is no ice!*/
     180                for(int i=0;i<ngrids_ocean;i++) if(icemask_oceangrid[i]>0.) icebase_oceangrid[i]=+9999.;
     181                xDelete<IssmDouble>(icemask_oceangrid);
     182                       
    160183                if(my_rank==0){
    161184                        ISSM_MPI_Send(icebase_oceangrid,ngrids_ocean,ISSM_MPI_DOUBLE,0,10001008,tomitgcmcomm);
     
    170193                xDelete<IssmDouble>(y_ice);
    171194                xDelete<IssmDouble>(z_ice);
    172                 xDelete<IssmDouble>(icebase);
    173195                xDelete<IssmDouble>(icebase_oceangrid);
    174196                xDelete<IssmDouble>(oceangridx);
     
    266288                        /*Interpolate ice base and mask onto ocean grid*/
    267289                        GetVectorFromInputsx(&icebase,femmodel,BaseEnum,VertexSIdEnum);
     290                        Options* options = new Options();
     291                        GenericOption<double> *odouble = new GenericOption<double>();
     292                        const char* name = "default";
     293                        odouble->name =xNew<char>(strlen(name)+1);
     294                        memcpy(odouble->name,name,(strlen(name)+1)*sizeof(char));
     295                        odouble->value=+9999.;
     296                        odouble->numel=1;
     297                        odouble->ndims=1;
     298                        odouble->size=NULL;
     299                        options->AddOption(odouble);
    268300                        InterpFromMeshToMesh2dx(&icebase_oceangrid,index_ice,lon_ice,lat_ice,ngrids_ice,nels_ice,
    269                                                 icebase,ngrids_ice,1,oceangridx,oceangridy,ngrids_ocean,NULL);
     301                                                icebase,ngrids_ice,1,oceangridx,oceangridy,ngrids_ocean,options);
     302                        delete options;
    270303                        xDelete<IssmDouble>(icebase);
    271304
    272305                        GetVectorFromInputsx(&icemask,femmodel,MaskIceLevelsetEnum,VertexSIdEnum);
     306                        Options* options2 = new Options();
     307                        GenericOption<double> *odouble2 = new GenericOption<double>();
     308                        const char* name2 = "default";
     309                        odouble2->name =xNew<char>(strlen(name2)+1);
     310                        memcpy(odouble2->name,name2,(strlen(name2)+1)*sizeof(char));
     311                        odouble2->value=+1.;
     312                        odouble2->numel=1;
     313                        odouble2->ndims=1;
     314                        odouble2->size=NULL;
     315                        options2->AddOption(odouble2);
    273316                        InterpFromMeshToMesh2dx(&icemask_oceangrid,index_ice,lon_ice,lat_ice,ngrids_ice,nels_ice,
    274                                                 icemask,ngrids_ice,1,oceangridx,oceangridy,ngrids_ocean,NULL);
     317                                                icemask,ngrids_ice,1,oceangridx,oceangridy,ngrids_ocean,options2);
     318                        delete options2;
    275319                        xDelete<IssmDouble>(icemask);
    276320
Note: See TracChangeset for help on using the changeset viewer.