- Timestamp:
- 03/09/18 08:44:31 (7 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk-jpl/src/c/modules/ModelProcessorx/CreateOutputDefinitions.cpp
r22410 r22507 65 65 else if (output_definition_enums[i]==MisfitEnum){ 66 66 /*Deal with misfits: {{{*/ 67 67 68 68 /*misfit variables: */ 69 69 int nummisfits; … … 157 157 /*}}}*/ 158 158 } 159 else if (output_definition_enums[i]==CfsurfacesquareEnum){ 160 /*Deal with cfsurfacesquare: {{{*/ 161 162 /*cfsurfacesquare variables: */ 163 int num_cfsurfacesquares; 164 char** cfsurfacesquare_name_s = NULL; 165 char** cfsurfacesquare_definitionstring_s = NULL; 166 char** cfsurfacesquare_model_string_s = NULL; 167 IssmDouble** cfsurfacesquare_observation_s = NULL; 168 char** cfsurfacesquare_observation_string_s = NULL; 169 int* cfsurfacesquare_observation_M_s = NULL; 170 int* cfsurfacesquare_observation_N_s = NULL; 171 IssmDouble** cfsurfacesquare_weights_s = NULL; 172 int* cfsurfacesquare_weights_M_s = NULL; 173 int* cfsurfacesquare_weights_N_s = NULL; 174 char** cfsurfacesquare_weights_string_s = NULL; 175 int* cfsurfacesquare_datatime_s = NULL; 176 177 /*Fetch name, model_string, observation, observation_string, etc ... (see src/m/classes/cfsurfacesquare.m): */ 178 iomodel->FetchMultipleData(&cfsurfacesquare_name_s,&num_cfsurfacesquares, "md.cfsurfacesquare.name"); 179 iomodel->FetchMultipleData(&cfsurfacesquare_definitionstring_s,&num_cfsurfacesquares, "md.cfsurfacesquare.definitionstring"); 180 iomodel->FetchMultipleData(&cfsurfacesquare_model_string_s,&num_cfsurfacesquares, "md.cfsurfacesquare.model_string"); 181 iomodel->FetchMultipleData(&cfsurfacesquare_observation_s,&cfsurfacesquare_observation_M_s,&cfsurfacesquare_observation_N_s,&num_cfsurfacesquares, "md.cfsurfacesquare.observation"); 182 iomodel->FetchMultipleData(&cfsurfacesquare_observation_string_s,&num_cfsurfacesquares, "md.cfsurfacesquare.observation_string"); 183 iomodel->FetchMultipleData(&cfsurfacesquare_weights_s,&cfsurfacesquare_weights_M_s,&cfsurfacesquare_weights_N_s,&num_cfsurfacesquares, "md.cfsurfacesquare.weights"); 184 iomodel->FetchMultipleData(&cfsurfacesquare_weights_string_s,&num_cfsurfacesquares, "md.cfsurfacesquare.weights_string"); 185 iomodel->FetchMultipleData(&cfsurfacesquare_datatime_s,&num_cfsurfacesquares, "md.cfsurfacesquare.datatime"); 186 187 for(j=0;j<num_cfsurfacesquares;j++){ 188 189 int obs_vector_type=0; 190 if ((cfsurfacesquare_observation_M_s[j]==iomodel->numberofvertices) || (cfsurfacesquare_observation_M_s[j]==iomodel->numberofvertices+1)){ 191 obs_vector_type=1; 192 } 193 else if ((cfsurfacesquare_observation_M_s[j]==iomodel->numberofelements) || (cfsurfacesquare_observation_M_s[j]==iomodel->numberofelements+1)){ 194 obs_vector_type=2; 195 } 196 else 197 _error_("cfsurfacesquare observation size not supported yet"); 198 199 int weight_vector_type=0; 200 if ((cfsurfacesquare_weights_M_s[j]==iomodel->numberofvertices) || (cfsurfacesquare_weights_M_s[j]==iomodel->numberofvertices+1)){ 201 weight_vector_type=1; 202 } 203 else if ((cfsurfacesquare_weights_M_s[j]==iomodel->numberofelements) || (cfsurfacesquare_weights_M_s[j]==iomodel->numberofelements+1)){ 204 weight_vector_type=2; 205 } 206 else 207 _error_("cfsurfacesquare weight size not supported yet"); 208 209 /*First create a cfsurfacesquare object for that specific string (cfsurfacesquare_model_string_s[j]):*/ 210 output_definitions->AddObject(new Cfsurfacesquare(cfsurfacesquare_name_s[j],StringToEnumx(cfsurfacesquare_definitionstring_s[j]),StringToEnumx(cfsurfacesquare_model_string_s[j]),StringToEnumx(cfsurfacesquare_observation_string_s[j]),StringToEnumx(cfsurfacesquare_weights_string_s[j]),cfsurfacesquare_datatime_s[j],false)); 211 212 /*Now, for this particular cfsurfacesquare object, make sure we plug into the elements: the observation, and the weights.*/ 213 for(int k=0;k<elements->Size();k++){ 214 215 Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(k)); 216 217 element->DatasetInputAdd(StringToEnumx(cfsurfacesquare_definitionstring_s[j]),cfsurfacesquare_observation_s[j], iomodel,cfsurfacesquare_observation_M_s[j],cfsurfacesquare_observation_N_s[j],obs_vector_type,StringToEnumx(cfsurfacesquare_observation_string_s[j]),7,SurfaceObservationEnum); 218 element->DatasetInputAdd(StringToEnumx(cfsurfacesquare_definitionstring_s[j]),cfsurfacesquare_weights_s[j], iomodel,cfsurfacesquare_weights_M_s[j],cfsurfacesquare_weights_N_s[j],weight_vector_type,StringToEnumx(cfsurfacesquare_weights_string_s[j]),7,WeightsSurfaceObservationEnum); 219 220 } 221 222 } 223 224 /*Free ressources:*/ 225 for(j=0;j<num_cfsurfacesquares;j++){ 226 char* string=NULL; 227 IssmDouble* matrix = NULL; 228 229 string = cfsurfacesquare_definitionstring_s[j]; xDelete<char>(string); 230 string = cfsurfacesquare_observation_string_s[j]; xDelete<char>(string); 231 string = cfsurfacesquare_model_string_s[j]; xDelete<char>(string); 232 string = cfsurfacesquare_weights_string_s[j]; xDelete<char>(string); 233 string = cfsurfacesquare_name_s[j]; xDelete<char>(string); 234 matrix = cfsurfacesquare_observation_s[j]; xDelete<IssmDouble>(matrix); 235 matrix = cfsurfacesquare_weights_s[j]; xDelete<IssmDouble>(matrix); 236 } 237 xDelete<char*>(cfsurfacesquare_name_s); 238 xDelete<char*>(cfsurfacesquare_model_string_s); 239 xDelete<char*>(cfsurfacesquare_definitionstring_s); 240 xDelete<IssmDouble*>(cfsurfacesquare_observation_s); 241 xDelete<char*>(cfsurfacesquare_observation_string_s); 242 xDelete<int>(cfsurfacesquare_observation_M_s); 243 xDelete<int>(cfsurfacesquare_observation_N_s); 244 xDelete<IssmDouble*>(cfsurfacesquare_weights_s); 245 xDelete<int>(cfsurfacesquare_weights_M_s); 246 xDelete<int>(cfsurfacesquare_weights_N_s); 247 xDelete<char*>(cfsurfacesquare_weights_string_s); 248 xDelete<int>(cfsurfacesquare_datatime_s); 249 /*}}}*/ 250 } 251 else if (output_definition_enums[i]==CfsurfacelogvelEnum){ 252 /*Deal with cfsurfacelogvel: {{{*/ 253 254 /*cfsurfacelogvel variables: */ 255 int num_cfsurfacelogvels; 256 char** cfsurfacelogvel_name = NULL; 257 char** cfsurfacelogvel_definitionstring = NULL; 258 IssmDouble** cfsurfacelogvel_vxobs = NULL; 259 IssmDouble** cfsurfacelogvel_vyobs = NULL; 260 char** cfsurfacelogvel_vxobs_string = NULL; 261 char** cfsurfacelogvel_vyobs_string = NULL; 262 int* cfsurfacelogvel_observation_M = NULL; 263 int* cfsurfacelogvel_observation_N = NULL; 264 IssmDouble** cfsurfacelogvel_weights = NULL; 265 int* cfsurfacelogvel_weights_M = NULL; 266 int* cfsurfacelogvel_weights_N = NULL; 267 char** cfsurfacelogvel_weightstring = NULL; 268 int* cfsurfacelogvel_datatime = NULL; 269 270 /*Fetch name, modeltring, observation, observationtring, etc ... (see src/m/classes/cfsurfacelogvel.m): */ 271 iomodel->FetchMultipleData(&cfsurfacelogvel_name,&num_cfsurfacelogvels, "md.cfsurfacelogvel.name"); 272 iomodel->FetchMultipleData(&cfsurfacelogvel_definitionstring,&num_cfsurfacelogvels, "md.cfsurfacelogvel.definitionstring"); 273 iomodel->FetchMultipleData(&cfsurfacelogvel_vxobs,&cfsurfacelogvel_observation_M,&cfsurfacelogvel_observation_N,&num_cfsurfacelogvels, "md.cfsurfacelogvel.vxobs"); 274 iomodel->FetchMultipleData(&cfsurfacelogvel_vxobs_string,&num_cfsurfacelogvels, "md.cfsurfacelogvel.vxobs_string"); 275 iomodel->FetchMultipleData(&cfsurfacelogvel_vyobs,&cfsurfacelogvel_observation_M,&cfsurfacelogvel_observation_N,&num_cfsurfacelogvels, "md.cfsurfacelogvel.vyobs"); 276 iomodel->FetchMultipleData(&cfsurfacelogvel_vyobs_string,&num_cfsurfacelogvels, "md.cfsurfacelogvel.vyobs_string"); iomodel->FetchMultipleData(&cfsurfacelogvel_weights,&cfsurfacelogvel_weights_M,&cfsurfacelogvel_weights_N,&num_cfsurfacelogvels, "md.cfsurfacelogvel.weights"); 277 iomodel->FetchMultipleData(&cfsurfacelogvel_weightstring,&num_cfsurfacelogvels, "md.cfsurfacelogvel.weights_string"); 278 _printf_("Num with weight string: "<<num_cfsurfacelogvels<<"\n"); 279 iomodel->FetchMultipleData(&cfsurfacelogvel_datatime,&num_cfsurfacelogvels, "md.cfsurfacelogvel.datatime"); 280 281 for(j=0;j<num_cfsurfacelogvels;j++){ 282 283 int obs_vector_type=0; 284 if ((cfsurfacelogvel_observation_M[j]==iomodel->numberofvertices) || (cfsurfacelogvel_observation_M[j]==iomodel->numberofvertices+1)){ 285 obs_vector_type=1; 286 } 287 else if ((cfsurfacelogvel_observation_M[j]==iomodel->numberofelements) || (cfsurfacelogvel_observation_M[j]==iomodel->numberofelements+1)){ 288 obs_vector_type=2; 289 } 290 else 291 _error_("cfsurfacelogvel observation size not supported yet"); 292 293 int weight_vector_type=0; 294 if ((cfsurfacelogvel_weights_M[j]==iomodel->numberofvertices) || (cfsurfacelogvel_weights_M[j]==iomodel->numberofvertices+1)){ 295 weight_vector_type=1; 296 } 297 else if ((cfsurfacelogvel_weights_M[j]==iomodel->numberofelements) || (cfsurfacelogvel_weights_M[j]==iomodel->numberofelements+1)){ 298 weight_vector_type=2; 299 } 300 else 301 _error_("cfsurfacelogvel weight size not supported yet"); 302 303 /*First create a cfsurfacelogvel object for that specific string (cfsurfacelogvel_modeltring[j]):*/ 304 output_definitions->AddObject(new Cfsurfacelogvel(cfsurfacelogvel_name[j],StringToEnumx(cfsurfacelogvel_definitionstring[j]),cfsurfacelogvel_datatime[j],false)); 305 306 /*Now, for this particular cfsurfacelogvel object, make sure we plug into the elements: the observation, and the weights.*/ 307 for(int k=0;k<elements->Size();k++){ 308 309 Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(k)); 310 311 element->DatasetInputAdd(StringToEnumx(cfsurfacelogvel_definitionstring[j]),cfsurfacelogvel_vxobs[j], iomodel,cfsurfacelogvel_observation_M[j],cfsurfacelogvel_observation_N[j],obs_vector_type,StringToEnumx(cfsurfacelogvel_vxobs_string[j]),7,VxObsEnum); 312 element->DatasetInputAdd(StringToEnumx(cfsurfacelogvel_definitionstring[j]),cfsurfacelogvel_vyobs[j], iomodel,cfsurfacelogvel_observation_M[j],cfsurfacelogvel_observation_N[j],obs_vector_type,StringToEnumx(cfsurfacelogvel_vyobs_string[j]),7,VyObsEnum); 313 element->DatasetInputAdd(StringToEnumx(cfsurfacelogvel_definitionstring[j]),cfsurfacelogvel_weights[j], iomodel,cfsurfacelogvel_weights_M[j],cfsurfacelogvel_weights_N[j],weight_vector_type,StringToEnumx(cfsurfacelogvel_weightstring[j]),7,WeightsSurfaceObservationEnum); 314 315 } 316 317 } 318 319 /*Free ressources:*/ 320 for(j=0;j<num_cfsurfacelogvels;j++){ 321 char* string=NULL; 322 IssmDouble* matrix = NULL; 323 324 string = cfsurfacelogvel_definitionstring[j]; xDelete<char>(string); 325 string = cfsurfacelogvel_vxobs_string[j]; xDelete<char>(string); 326 string = cfsurfacelogvel_vyobs_string[j]; xDelete<char>(string); 327 string = cfsurfacelogvel_weightstring[j]; xDelete<char>(string); 328 string = cfsurfacelogvel_name[j]; xDelete<char>(string); 329 matrix = cfsurfacelogvel_weights[j]; xDelete<IssmDouble>(matrix); 330 matrix = cfsurfacelogvel_vxobs[j]; xDelete<IssmDouble>(matrix); 331 matrix = cfsurfacelogvel_vyobs[j]; xDelete<IssmDouble>(matrix); 332 } 333 xDelete<char*>(cfsurfacelogvel_name); 334 xDelete<char*>(cfsurfacelogvel_definitionstring); 335 xDelete<int>(cfsurfacelogvel_observation_M); 336 xDelete<IssmDouble*>(cfsurfacelogvel_vxobs); 337 xDelete<IssmDouble*>(cfsurfacelogvel_vyobs); 338 xDelete<char*>(cfsurfacelogvel_vxobs_string); 339 xDelete<char*>(cfsurfacelogvel_vyobs_string); 340 xDelete<int>(cfsurfacelogvel_observation_N); 341 xDelete<IssmDouble*>(cfsurfacelogvel_weights); 342 xDelete<int>(cfsurfacelogvel_weights_M); 343 xDelete<int>(cfsurfacelogvel_weights_N); 344 xDelete<char*>(cfsurfacelogvel_weightstring); 345 xDelete<int>(cfsurfacelogvel_datatime); 346 /*}}}*/ 347 } 159 348 else if (output_definition_enums[i]==NodalvalueEnum){ 160 349 /*Deal with nodal values: {{{*/
Note:
See TracChangeset
for help on using the changeset viewer.