Changeset 27814
- Timestamp:
- 06/30/23 14:01:03 (21 months ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
TabularUnified issm/trunk-jpl/src/c/modules/OceanExchangeDatax/OceanExchangeDatax.cpp ¶
r27812 r27814 28 28 IssmDouble *oceangridx; 29 29 IssmDouble *oceangridy; 30 IssmDouble *ice base_oceangrid = NULL;30 IssmDouble *icethickness_oceangrid = NULL; 31 31 IssmDouble *icemask_oceangrid = NULL; 32 32 IssmDouble* x_ice = NULL; … … 34 34 IssmDouble* lat_ice = NULL; 35 35 IssmDouble* lon_ice = NULL; 36 IssmDouble* ice base= NULL;36 IssmDouble* icethickness = NULL; 37 37 IssmDouble* icemask = NULL; 38 38 IssmDouble* melt_mesh = NULL; … … 81 81 } 82 82 83 /*Interpolate ice baseand mask onto ocean grid*/83 /*Interpolate ice thickness and mask onto ocean grid*/ 84 84 femmodel->GetMesh(femmodel->vertices,femmodel->elements,&x_ice,&y_ice,&index_ice); 85 85 BamgTriangulatex(&index_ocean,&nels_ocean,oceangridx,oceangridy,ngrids_ocean); 86 86 femmodel->vertices->XYList(&lon_ice,&lat_ice); 87 87 //femmodel->vertices->LatLonList(&lat_ice,&lon_ice); 88 GetVectorFromInputsx(&ice base,femmodel,BaseEnum,VertexSIdEnum);88 GetVectorFromInputsx(&icethickness,femmodel,ThicknessEnum,VertexSIdEnum); 89 89 Options* options = new Options(); 90 90 GenericOption<double> *odouble = new GenericOption<double>(); … … 96 96 odouble->size[1]=1; 97 97 options->AddOption(odouble); 98 InterpFromMeshToMesh2dx(&ice base_oceangrid,index_ice,lon_ice,lat_ice,ngrids_ice,nels_ice,99 ice base,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); 100 100 delete options; 101 xDelete<IssmDouble>(ice base);101 xDelete<IssmDouble>(icethickness); 102 102 103 103 GetVectorFromInputsx(&icemask,femmodel,MaskIceLevelsetEnum,VertexSIdEnum); … … 117 117 118 118 /*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.; 120 122 xDelete<IssmDouble>(icemask_oceangrid); 121 123 122 if(init_stage==true){ //just send ice base124 if(init_stage==true){ //just send icethickness 123 125 if(my_rank==0){ 124 ISSM_MPI_Send(ice base_oceangrid,ngrids_ocean,ISSM_MPI_DOUBLE,0,10001008,tomitgcmcomm);126 ISSM_MPI_Send(icethickness_oceangrid,ngrids_ocean,ISSM_MPI_DOUBLE,0,10001008,tomitgcmcomm); 125 127 } 126 128 } 127 129 else{ //send and receive exchanged data 128 femmodel->parameters->FindParam(&rho_ice,MaterialsRhoIceEnum);129 130 femmodel->parameters->FindParam(&yts,ConstantsYtsEnum); 130 131 if(my_rank==0){ … … 134 135 oceanmelt = xNew<IssmDouble>(ngrids_ocean); 135 136 ISSM_MPI_Recv(oceanmelt,ngrids_ocean,ISSM_MPI_DOUBLE,0,10001007,tomitgcmcomm,&status); 136 ISSM_MPI_Send(ice base_oceangrid,ngrids_ocean,ISSM_MPI_DOUBLE,0,10001008,tomitgcmcomm);137 ISSM_MPI_Send(icethickness_oceangrid,ngrids_ocean,ISSM_MPI_DOUBLE,0,10001008,tomitgcmcomm); 137 138 } 138 139 ISSM_MPI_Bcast(&oceantime,1,ISSM_MPI_DOUBLE,0,IssmComm::GetComm()); … … 156 157 xDelete<IssmDouble>(x_ice); 157 158 xDelete<IssmDouble>(y_ice); 158 xDelete<IssmDouble>(ice base_oceangrid);159 xDelete<IssmDouble>(icethickness_oceangrid); 159 160 xDelete<IssmDouble>(oceangridx); 160 161 xDelete<IssmDouble>(oceangridy);
Note:
See TracChangeset
for help on using the changeset viewer.