Ice Sheet System Model  4.18
Code documentation
Functions
CreateOutputDefinitions.cpp File Reference
#include "../../toolkits/toolkits.h"
#include "../../classes/classes.h"
#include "../../shared/shared.h"
#include "./ModelProcessorx.h"

Go to the source code of this file.

Functions

void CreateOutputDefinitions (Elements *elements, Parameters *parameters, Inputs2 *inputs2, IoModel *iomodel)
 

Function Documentation

◆ CreateOutputDefinitions()

void CreateOutputDefinitions ( Elements elements,
Parameters parameters,
Inputs2 inputs2,
IoModel iomodel 
)

Definition at line 10 of file CreateOutputDefinitions.cpp.

10  {
11 
12  int i,j;
13 
14  DataSet* output_definitions = NULL;
15  int* output_definition_enums = NULL;
16  int num_output_definitions;
17 
18  /*Create output_definitions dataset: */
19  output_definitions=new DataSet();
20  char** out_strings = NULL;
21  iomodel->FetchData(&out_strings,&num_output_definitions,"md.outputdefinition.list");
22  if(num_output_definitions>0){
23  output_definition_enums=xNew<int>(num_output_definitions);
24  for(int i=0;i<num_output_definitions;i++){
25  output_definition_enums[i]=StringToEnumx(out_strings[i]);
26  }
27  }
28  // free data:
29  for(int i=0;i<num_output_definitions;i++) xDelete<char>(out_strings[i]);
30  xDelete<char*>(out_strings);
31 
32  if(num_output_definitions){
33  for (i=0;i<num_output_definitions;i++){
34  if (output_definition_enums[i]==MassfluxatgateEnum){
35  /*Deal with mass flux gates:{{{ */
36 
37  /*massfluxatgate variables: */
38  int temp,numgates;
39  char **gatenames = NULL;
40  char **gatedefinitionstrings = NULL;
41  IssmDouble **gatesegments = NULL;
42  int *gatesegments_M = NULL;
43 
44  /*Fetch segments and names: */
45  iomodel->FetchMultipleData(&gatenames,&numgates, "md.massfluxatgate.name");
46  iomodel->FetchMultipleData(&gatedefinitionstrings,&temp, "md.massfluxatgate.definitionstring"); _assert_(temp==numgates);
47  iomodel->FetchMultipleData(&gatesegments,&gatesegments_M,NULL,&temp, "md.massfluxatgate.segments"); _assert_(temp==numgates);
48 
49  for(j=0;j<numgates;j++){
50  output_definitions->AddObject(new Massfluxatgate<IssmDouble>(gatenames[j],StringToEnumx(gatedefinitionstrings[j]),gatesegments_M[j],gatesegments[j]));
51  }
52  /*Free ressources:*/
53  for(j=0;j<numgates;j++){
54  char* string = gatenames[j]; xDelete<char>(string);
55  char* string2 = gatedefinitionstrings[j]; xDelete<char>(string2);
56  IssmDouble* gate = gatesegments[j]; xDelete<IssmDouble>(gate);
57  }
58  xDelete<char*>(gatenames);
59  xDelete<IssmDouble*>(gatesegments);
60  xDelete<int>(gatesegments_M);
61  xDelete<char*>(gatedefinitionstrings);
62  /*}}}*/
63  }
64  else if (output_definition_enums[i]==MisfitEnum){
65  /*Deal with misfits: {{{*/
66 
67  /*misfit variables: */
68  int nummisfits;
69  char** misfit_name_s = NULL;
70  char** misfit_definitionstring_s = NULL;
71  char** misfit_model_string_s = NULL;
72  IssmDouble** misfit_observation_s = NULL;
73  char** misfit_observation_string_s = NULL;
74  int* misfit_observation_M_s = NULL;
75  int* misfit_observation_N_s = NULL;
76  int* misfit_local_s = NULL;
77  char** misfit_timeinterpolation_s = NULL;
78  IssmDouble** misfit_weights_s = NULL;
79  int* misfit_weights_M_s = NULL;
80  int* misfit_weights_N_s = NULL;
81  char** misfit_weights_string_s = NULL;
82 
83  /*Fetch name, model_string, observation, observation_string, etc ... (see src/m/classes/misfit.m): */
84  iomodel->FetchMultipleData(&misfit_name_s,&nummisfits, "md.misfit.name");
85  iomodel->FetchMultipleData(&misfit_definitionstring_s,&nummisfits, "md.misfit.definitionstring");
86  iomodel->FetchMultipleData(&misfit_model_string_s,&nummisfits, "md.misfit.model_string");
87  iomodel->FetchMultipleData(&misfit_observation_s,&misfit_observation_M_s,&misfit_observation_N_s,&nummisfits, "md.misfit.observation");
88  iomodel->FetchMultipleData(&misfit_observation_string_s,&nummisfits, "md.misfit.observation_string");
89  iomodel->FetchMultipleData(&misfit_timeinterpolation_s,&nummisfits, "md.misfit.timeinterpolation");
90  iomodel->FetchMultipleData(&misfit_local_s,&nummisfits, "md.misfit.local");
91  iomodel->FetchMultipleData(&misfit_weights_s,&misfit_weights_M_s,&misfit_weights_N_s,&nummisfits, "md.misfit.weights");
92  iomodel->FetchMultipleData(&misfit_weights_string_s,&nummisfits, "md.misfit.weights_string");
93 
94  for(j=0;j<nummisfits;j++){
95 
96  int obs_vector_type=0;
97  if ((misfit_observation_M_s[j]==iomodel->numberofvertices) || (misfit_observation_M_s[j]==iomodel->numberofvertices+1)){
98  obs_vector_type=1;
99  }
100  else if ((misfit_observation_M_s[j]==iomodel->numberofelements) || (misfit_observation_M_s[j]==iomodel->numberofelements+1)){
101  obs_vector_type=2;
102  }
103  else
104  _error_("misfit observation size not supported yet");
105 
106  int weight_vector_type=0;
107  if ((misfit_weights_M_s[j]==iomodel->numberofvertices) || (misfit_weights_M_s[j]==iomodel->numberofvertices+1)){
108  weight_vector_type=1;
109  }
110  else if ((misfit_weights_M_s[j]==iomodel->numberofelements) || (misfit_weights_M_s[j]==iomodel->numberofelements+1)){
111  weight_vector_type=2;
112  }
113  else
114  _error_("misfit weight size not supported yet");
115 
116  /*First create a misfit object for that specific string (misfit_model_string_s[j]):*/
117  output_definitions->AddObject(new Misfit(misfit_name_s[j],StringToEnumx(misfit_definitionstring_s[j]),StringToEnumx(misfit_model_string_s[j]),StringToEnumx(misfit_observation_string_s[j]),misfit_timeinterpolation_s[j],misfit_local_s[j],StringToEnumx(misfit_weights_string_s[j])));
118 
119  /*Now, for this particular misfit object, make sure we plug into the elements: the observation, and the weights.*/
120  for(int k=0;k<elements->Size();k++){
121  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(k));
122  element->InputCreate(misfit_observation_s[j],inputs2,iomodel,misfit_observation_M_s[j],misfit_observation_N_s[j],obs_vector_type,StringToEnumx(misfit_observation_string_s[j]),7);
123  element->InputCreate(misfit_weights_s[j],inputs2,iomodel,misfit_weights_M_s[j],misfit_weights_N_s[j],weight_vector_type,StringToEnumx(misfit_weights_string_s[j]),7);
124  }
125 
126  }
127 
128  /*Free ressources:*/
129  for(j=0;j<nummisfits;j++){
130  char* string=NULL;
131  IssmDouble* matrix = NULL;
132  string = misfit_definitionstring_s[j]; xDelete<char>(string);
133  string = misfit_observation_string_s[j]; xDelete<char>(string);
134  string = misfit_model_string_s[j]; xDelete<char>(string);
135  string = misfit_weights_string_s[j]; xDelete<char>(string);
136  string = misfit_name_s[j]; xDelete<char>(string);
137  string = misfit_timeinterpolation_s[j]; xDelete<char>(string);
138  matrix = misfit_observation_s[j]; xDelete<IssmDouble>(matrix);
139  matrix = misfit_weights_s[j]; xDelete<IssmDouble>(matrix);
140  }
141  xDelete<char*>(misfit_name_s);
142  xDelete<char*>(misfit_model_string_s);
143  xDelete<char*>(misfit_definitionstring_s);
144  xDelete<IssmDouble*>(misfit_observation_s);
145  xDelete<char*>(misfit_observation_string_s);
146  xDelete<int>(misfit_observation_M_s);
147  xDelete<int>(misfit_observation_N_s);
148  xDelete<int>(misfit_local_s);
149  xDelete<char*>(misfit_timeinterpolation_s);
150  xDelete<IssmDouble*>(misfit_weights_s);
151  xDelete<int>(misfit_weights_M_s);
152  xDelete<int>(misfit_weights_N_s);
153  xDelete<char*>(misfit_weights_string_s);
154  /*}}}*/
155  }
156  else if (output_definition_enums[i]==CfsurfacesquareEnum){
157  /*Deal with cfsurfacesquare: {{{*/
158 
159  /*cfsurfacesquare variables: */
160  int num_cfsurfacesquares;
161  char** cfsurfacesquare_name_s = NULL;
162  char** cfsurfacesquare_definitionstring_s = NULL;
163  char** cfsurfacesquare_model_string_s = NULL;
164  IssmDouble** cfsurfacesquare_observation_s = NULL;
165  char** cfsurfacesquare_observation_string_s = NULL;
166  int* cfsurfacesquare_observation_M_s = NULL;
167  int* cfsurfacesquare_observation_N_s = NULL;
168  IssmDouble** cfsurfacesquare_weights_s = NULL;
169  int* cfsurfacesquare_weights_M_s = NULL;
170  int* cfsurfacesquare_weights_N_s = NULL;
171  char** cfsurfacesquare_weights_string_s = NULL;
172  int* cfsurfacesquare_datatime_s = NULL;
173 
174  /*Fetch name, model_string, observation, observation_string, etc ... (see src/m/classes/cfsurfacesquare.m): */
175  iomodel->FetchMultipleData(&cfsurfacesquare_name_s,&num_cfsurfacesquares, "md.cfsurfacesquare.name");
176  iomodel->FetchMultipleData(&cfsurfacesquare_definitionstring_s,&num_cfsurfacesquares, "md.cfsurfacesquare.definitionstring");
177  iomodel->FetchMultipleData(&cfsurfacesquare_model_string_s,&num_cfsurfacesquares, "md.cfsurfacesquare.model_string");
178  iomodel->FetchMultipleData(&cfsurfacesquare_observation_s,&cfsurfacesquare_observation_M_s,&cfsurfacesquare_observation_N_s,&num_cfsurfacesquares, "md.cfsurfacesquare.observation");
179  iomodel->FetchMultipleData(&cfsurfacesquare_observation_string_s,&num_cfsurfacesquares, "md.cfsurfacesquare.observation_string");
180  iomodel->FetchMultipleData(&cfsurfacesquare_weights_s,&cfsurfacesquare_weights_M_s,&cfsurfacesquare_weights_N_s,&num_cfsurfacesquares, "md.cfsurfacesquare.weights");
181  iomodel->FetchMultipleData(&cfsurfacesquare_weights_string_s,&num_cfsurfacesquares, "md.cfsurfacesquare.weights_string");
182  iomodel->FetchMultipleData(&cfsurfacesquare_datatime_s,&num_cfsurfacesquares, "md.cfsurfacesquare.datatime");
183 
184  for(j=0;j<num_cfsurfacesquares;j++){
185 
186  int obs_vector_type=0;
187  if ((cfsurfacesquare_observation_M_s[j]==iomodel->numberofvertices) || (cfsurfacesquare_observation_M_s[j]==iomodel->numberofvertices+1)){
188  obs_vector_type=1;
189  }
190  else if ((cfsurfacesquare_observation_M_s[j]==iomodel->numberofelements) || (cfsurfacesquare_observation_M_s[j]==iomodel->numberofelements+1)){
191  obs_vector_type=2;
192  }
193  else
194  _error_("cfsurfacesquare observation size not supported yet");
195 
196  int weight_vector_type=0;
197  if ((cfsurfacesquare_weights_M_s[j]==iomodel->numberofvertices) || (cfsurfacesquare_weights_M_s[j]==iomodel->numberofvertices+1)){
198  weight_vector_type=1;
199  }
200  else if ((cfsurfacesquare_weights_M_s[j]==iomodel->numberofelements) || (cfsurfacesquare_weights_M_s[j]==iomodel->numberofelements+1)){
201  weight_vector_type=2;
202  }
203  else
204  _error_("cfsurfacesquare weight size not supported yet");
205 
206  /*First create a cfsurfacesquare object for that specific string (cfsurfacesquare_model_string_s[j]):*/
207  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));
208 
209  /*Now, for this particular cfsurfacesquare object, make sure we plug into the elements: the observation, and the weights.*/
210  for(int k=0;k<elements->Size();k++){
211  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(k));
212  element->DatasetInputAdd(StringToEnumx(cfsurfacesquare_definitionstring_s[j]),cfsurfacesquare_observation_s[j],inputs2,iomodel,cfsurfacesquare_observation_M_s[j],cfsurfacesquare_observation_N_s[j],obs_vector_type,StringToEnumx(cfsurfacesquare_observation_string_s[j]),7,SurfaceObservationEnum);
213  element->DatasetInputAdd(StringToEnumx(cfsurfacesquare_definitionstring_s[j]),cfsurfacesquare_weights_s[j],inputs2,iomodel,cfsurfacesquare_weights_M_s[j],cfsurfacesquare_weights_N_s[j],weight_vector_type,StringToEnumx(cfsurfacesquare_weights_string_s[j]),7,WeightsSurfaceObservationEnum);
214 
215  }
216 
217  }
218 
219  /*Free ressources:*/
220  for(j=0;j<num_cfsurfacesquares;j++){
221  char* string=NULL;
222  IssmDouble* matrix = NULL;
223 
224  string = cfsurfacesquare_definitionstring_s[j]; xDelete<char>(string);
225  string = cfsurfacesquare_observation_string_s[j]; xDelete<char>(string);
226  string = cfsurfacesquare_model_string_s[j]; xDelete<char>(string);
227  string = cfsurfacesquare_weights_string_s[j]; xDelete<char>(string);
228  string = cfsurfacesquare_name_s[j]; xDelete<char>(string);
229  matrix = cfsurfacesquare_observation_s[j]; xDelete<IssmDouble>(matrix);
230  matrix = cfsurfacesquare_weights_s[j]; xDelete<IssmDouble>(matrix);
231  }
232  xDelete<char*>(cfsurfacesquare_name_s);
233  xDelete<char*>(cfsurfacesquare_model_string_s);
234  xDelete<char*>(cfsurfacesquare_definitionstring_s);
235  xDelete<IssmDouble*>(cfsurfacesquare_observation_s);
236  xDelete<char*>(cfsurfacesquare_observation_string_s);
237  xDelete<int>(cfsurfacesquare_observation_M_s);
238  xDelete<int>(cfsurfacesquare_observation_N_s);
239  xDelete<IssmDouble*>(cfsurfacesquare_weights_s);
240  xDelete<int>(cfsurfacesquare_weights_M_s);
241  xDelete<int>(cfsurfacesquare_weights_N_s);
242  xDelete<char*>(cfsurfacesquare_weights_string_s);
243  xDelete<int>(cfsurfacesquare_datatime_s);
244  /*}}}*/
245  }
246  else if (output_definition_enums[i]==CfdragcoeffabsgradEnum){
247  /*Deal with cfdragcoeffabsgrad: {{{*/
248 
249  /*cfdragcoeffabsgrad variables: */
250  int num_cfdragcoeffabsgrads;
251  char** cfdragcoeffabsgrad_name_s = NULL;
252  char** cfdragcoeffabsgrad_definitionstring_s = NULL;
253  IssmDouble** cfdragcoeffabsgrad_weights_s = NULL;
254  int* cfdragcoeffabsgrad_weights_M_s = NULL;
255  int* cfdragcoeffabsgrad_weights_N_s = NULL;
256  char** cfdragcoeffabsgrad_weights_string_s = NULL;
257 
258  /*Fetch name, model_string, observation, observation_string, etc ... (see src/m/classes/cfdragcoeffabsgrad.m): */
259  iomodel->FetchMultipleData(&cfdragcoeffabsgrad_name_s,&num_cfdragcoeffabsgrads, "md.cfdragcoeffabsgrad.name");
260  iomodel->FetchMultipleData(&cfdragcoeffabsgrad_definitionstring_s,&num_cfdragcoeffabsgrads, "md.cfdragcoeffabsgrad.definitionstring");
261  iomodel->FetchMultipleData(&cfdragcoeffabsgrad_weights_s,&cfdragcoeffabsgrad_weights_M_s,&cfdragcoeffabsgrad_weights_N_s,&num_cfdragcoeffabsgrads, "md.cfdragcoeffabsgrad.weights");
262  iomodel->FetchMultipleData(&cfdragcoeffabsgrad_weights_string_s,&num_cfdragcoeffabsgrads, "md.cfdragcoeffabsgrad.weights_string");
263 
264  for(j=0;j<num_cfdragcoeffabsgrads;j++){
265 
266  int weight_vector_type=0;
267  if ((cfdragcoeffabsgrad_weights_M_s[j]==iomodel->numberofvertices) || (cfdragcoeffabsgrad_weights_M_s[j]==iomodel->numberofvertices+1)){
268  weight_vector_type=1;
269  }
270  else if ((cfdragcoeffabsgrad_weights_M_s[j]==iomodel->numberofelements) || (cfdragcoeffabsgrad_weights_M_s[j]==iomodel->numberofelements+1)){
271  weight_vector_type=2;
272  }
273  else
274  _error_("cfdragcoeffabsgrad weight size not supported yet");
275 
276  /*First create a cfdragcoeffabsgrad object for that specific string (cfdragcoeffabsgrad_model_string_s[j]):*/
277  output_definitions->AddObject(new Cfdragcoeffabsgrad(cfdragcoeffabsgrad_name_s[j],StringToEnumx(cfdragcoeffabsgrad_definitionstring_s[j]),StringToEnumx(cfdragcoeffabsgrad_weights_string_s[j]),false));
278 
279  /*Now, for this particular cfdragcoeffabsgrad object, make sure we plug into the elements: the observation, and the weights.*/
280  for(int k=0;k<elements->Size();k++){
281 
282  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(k));
283 
284  element->DatasetInputAdd(StringToEnumx(cfdragcoeffabsgrad_definitionstring_s[j]),cfdragcoeffabsgrad_weights_s[j],inputs2,iomodel,cfdragcoeffabsgrad_weights_M_s[j],cfdragcoeffabsgrad_weights_N_s[j],weight_vector_type,StringToEnumx(cfdragcoeffabsgrad_weights_string_s[j]),7,WeightsSurfaceObservationEnum);
285 
286  }
287 
288  }
289 
290  /*Free ressources:*/
291  for(j=0;j<num_cfdragcoeffabsgrads;j++){
292  char* string=NULL;
293  IssmDouble* matrix = NULL;
294 
295  string = cfdragcoeffabsgrad_definitionstring_s[j]; xDelete<char>(string);
296  string = cfdragcoeffabsgrad_weights_string_s[j]; xDelete<char>(string);
297  string = cfdragcoeffabsgrad_name_s[j]; xDelete<char>(string);
298  matrix = cfdragcoeffabsgrad_weights_s[j]; xDelete<IssmDouble>(matrix);
299  }
300  xDelete<char*>(cfdragcoeffabsgrad_name_s);
301  xDelete<char*>(cfdragcoeffabsgrad_definitionstring_s);
302  xDelete<IssmDouble*>(cfdragcoeffabsgrad_weights_s);
303  xDelete<int>(cfdragcoeffabsgrad_weights_M_s);
304  xDelete<int>(cfdragcoeffabsgrad_weights_N_s);
305  xDelete<char*>(cfdragcoeffabsgrad_weights_string_s);
306  /*}}}*/
307  }
308  else if (output_definition_enums[i]==CfsurfacelogvelEnum){
309  /*Deal with cfsurfacelogvel: {{{*/
310 
311  /*cfsurfacelogvel variables: */
312  int num_cfsurfacelogvels;
313  char** cfsurfacelogvel_name = NULL;
314  char** cfsurfacelogvel_definitionstring = NULL;
315  IssmDouble** cfsurfacelogvel_vxobs = NULL;
316  IssmDouble** cfsurfacelogvel_vyobs = NULL;
317  char** cfsurfacelogvel_vxobs_string = NULL;
318  char** cfsurfacelogvel_vyobs_string = NULL;
319  int* cfsurfacelogvel_observation_M = NULL;
320  int* cfsurfacelogvel_observation_N = NULL;
321  IssmDouble** cfsurfacelogvel_weights = NULL;
322  int* cfsurfacelogvel_weights_M = NULL;
323  int* cfsurfacelogvel_weights_N = NULL;
324  char** cfsurfacelogvel_weightstring = NULL;
325  int* cfsurfacelogvel_datatime = NULL;
326 
327  /*Fetch name, modeltring, observation, observationtring, etc ... (see src/m/classes/cfsurfacelogvel.m): */
328  iomodel->FetchMultipleData(&cfsurfacelogvel_name,&num_cfsurfacelogvels, "md.cfsurfacelogvel.name");
329  iomodel->FetchMultipleData(&cfsurfacelogvel_definitionstring,&num_cfsurfacelogvels, "md.cfsurfacelogvel.definitionstring");
330  iomodel->FetchMultipleData(&cfsurfacelogvel_vxobs,&cfsurfacelogvel_observation_M,&cfsurfacelogvel_observation_N,&num_cfsurfacelogvels, "md.cfsurfacelogvel.vxobs");
331  iomodel->FetchMultipleData(&cfsurfacelogvel_vxobs_string,&num_cfsurfacelogvels, "md.cfsurfacelogvel.vxobs_string");
332  iomodel->FetchMultipleData(&cfsurfacelogvel_vyobs,&cfsurfacelogvel_observation_M,&cfsurfacelogvel_observation_N,&num_cfsurfacelogvels, "md.cfsurfacelogvel.vyobs");
333  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");
334  iomodel->FetchMultipleData(&cfsurfacelogvel_weightstring,&num_cfsurfacelogvels, "md.cfsurfacelogvel.weights_string");
335  iomodel->FetchMultipleData(&cfsurfacelogvel_datatime,&num_cfsurfacelogvels, "md.cfsurfacelogvel.datatime");
336 
337  for(j=0;j<num_cfsurfacelogvels;j++){
338 
339  int obs_vector_type=0;
340  if ((cfsurfacelogvel_observation_M[j]==iomodel->numberofvertices) || (cfsurfacelogvel_observation_M[j]==iomodel->numberofvertices+1)){
341  obs_vector_type=1;
342  }
343  else if ((cfsurfacelogvel_observation_M[j]==iomodel->numberofelements) || (cfsurfacelogvel_observation_M[j]==iomodel->numberofelements+1)){
344  obs_vector_type=2;
345  }
346  else
347  _error_("cfsurfacelogvel observation size not supported yet");
348 
349  int weight_vector_type=0;
350  if ((cfsurfacelogvel_weights_M[j]==iomodel->numberofvertices) || (cfsurfacelogvel_weights_M[j]==iomodel->numberofvertices+1)){
351  weight_vector_type=1;
352  }
353  else if ((cfsurfacelogvel_weights_M[j]==iomodel->numberofelements) || (cfsurfacelogvel_weights_M[j]==iomodel->numberofelements+1)){
354  weight_vector_type=2;
355  }
356  else
357  _error_("cfsurfacelogvel weight size not supported yet");
358 
359  /*First create a cfsurfacelogvel object for that specific string (cfsurfacelogvel_modeltring[j]):*/
360  output_definitions->AddObject(new Cfsurfacelogvel(cfsurfacelogvel_name[j],StringToEnumx(cfsurfacelogvel_definitionstring[j]),cfsurfacelogvel_datatime[j],false));
361 
362  /*Now, for this particular cfsurfacelogvel object, make sure we plug into the elements: the observation, and the weights.*/
363  for(int k=0;k<elements->Size();k++){
364 
365  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(k));
366 
367  element->DatasetInputAdd(StringToEnumx(cfsurfacelogvel_definitionstring[j]),cfsurfacelogvel_vxobs[j],inputs2,iomodel,cfsurfacelogvel_observation_M[j],cfsurfacelogvel_observation_N[j],obs_vector_type,StringToEnumx(cfsurfacelogvel_vxobs_string[j]),7,VxObsEnum);
368  element->DatasetInputAdd(StringToEnumx(cfsurfacelogvel_definitionstring[j]),cfsurfacelogvel_vyobs[j],inputs2,iomodel,cfsurfacelogvel_observation_M[j],cfsurfacelogvel_observation_N[j],obs_vector_type,StringToEnumx(cfsurfacelogvel_vyobs_string[j]),7,VyObsEnum);
369  element->DatasetInputAdd(StringToEnumx(cfsurfacelogvel_definitionstring[j]),cfsurfacelogvel_weights[j],inputs2,iomodel,cfsurfacelogvel_weights_M[j],cfsurfacelogvel_weights_N[j],weight_vector_type,StringToEnumx(cfsurfacelogvel_weightstring[j]),7,WeightsSurfaceObservationEnum);
370 
371  }
372 
373  }
374 
375  /*Free ressources:*/
376  for(j=0;j<num_cfsurfacelogvels;j++){
377  char* string=NULL;
378  IssmDouble* matrix = NULL;
379 
380  string = cfsurfacelogvel_definitionstring[j]; xDelete<char>(string);
381  string = cfsurfacelogvel_vxobs_string[j]; xDelete<char>(string);
382  string = cfsurfacelogvel_vyobs_string[j]; xDelete<char>(string);
383  string = cfsurfacelogvel_weightstring[j]; xDelete<char>(string);
384  string = cfsurfacelogvel_name[j]; xDelete<char>(string);
385  matrix = cfsurfacelogvel_weights[j]; xDelete<IssmDouble>(matrix);
386  matrix = cfsurfacelogvel_vxobs[j]; xDelete<IssmDouble>(matrix);
387  matrix = cfsurfacelogvel_vyobs[j]; xDelete<IssmDouble>(matrix);
388  }
389  xDelete<char*>(cfsurfacelogvel_name);
390  xDelete<char*>(cfsurfacelogvel_definitionstring);
391  xDelete<int>(cfsurfacelogvel_observation_M);
392  xDelete<IssmDouble*>(cfsurfacelogvel_vxobs);
393  xDelete<IssmDouble*>(cfsurfacelogvel_vyobs);
394  xDelete<char*>(cfsurfacelogvel_vxobs_string);
395  xDelete<char*>(cfsurfacelogvel_vyobs_string);
396  xDelete<int>(cfsurfacelogvel_observation_N);
397  xDelete<IssmDouble*>(cfsurfacelogvel_weights);
398  xDelete<int>(cfsurfacelogvel_weights_M);
399  xDelete<int>(cfsurfacelogvel_weights_N);
400  xDelete<char*>(cfsurfacelogvel_weightstring);
401  xDelete<int>(cfsurfacelogvel_datatime);
402  /*}}}*/
403  }
404  else if (output_definition_enums[i]==NodalvalueEnum){
405  /*Deal with nodal values: {{{*/
406 
407  /*nodal value variables: */
408  int numnodalvalues;
409  char** nodalvalue_name_s = NULL;
410  char** nodalvalue_definitionstrings = NULL;
411  char** nodalvalue_modelstrings = NULL;
412  int* nodalvalue_node_s = NULL;
413 
414  /*Fetch name, model_enum, etc ... (see src/m/classes/nodalvalue.m): */
415  iomodel->FetchMultipleData(&nodalvalue_name_s,&numnodalvalues, "md.nodalvalue.name");
416  iomodel->FetchMultipleData(&nodalvalue_definitionstrings,&numnodalvalues, "md.nodalvalue.definitionenum");
417  iomodel->FetchMultipleData(&nodalvalue_modelstrings,&numnodalvalues, "md.nodalvalue.model_enum");
418  iomodel->FetchMultipleData(&nodalvalue_node_s,&numnodalvalues, "md.nodalvalue.node");
419 
420  for(j=0;j<numnodalvalues;j++){
421 
422  /*First create a nodalvalue object for that specific enum (nodalvalue_model_enum_s[j]):*/
423  output_definitions->AddObject(new Nodalvalue(nodalvalue_name_s[j],StringToEnumx(nodalvalue_definitionstrings[j]),StringToEnumx(nodalvalue_modelstrings[j]),nodalvalue_node_s[j]-1)); //-1 because matlab to c indexing.
424  }
425 
426  /*Free ressources:*/
427  for(j=0;j<numnodalvalues;j++){
428  char* string=NULL;
429  string = nodalvalue_name_s[j]; xDelete<char>(string);
430  }
431  xDelete<char*>(nodalvalue_name_s);
432  xDelete<char*>(nodalvalue_modelstrings);
433  xDelete<char*>(nodalvalue_definitionstrings);
434  xDelete<int>(nodalvalue_node_s);
435  /*}}}*/
436  }
437  else if (output_definition_enums[i]==MassconEnum){
438  /*Deal with masscons: {{{*/
439 
440  /*masscon variables: */
441  int nummasscons;
442  char** masscon_name_s = NULL;
443  char** masscon_definitionstring_s = NULL;
444  IssmDouble** masscon_levelset_s = NULL;
445  int* masscon_levelset_M_s = NULL;
446  int* masscon_levelset_N_s = NULL;
447 
448  /*Fetch name and levelset, etc ... (see src/m/classes/masscon.m): */
449  iomodel->FetchMultipleData(&masscon_name_s,&nummasscons, "md.masscon.name");
450  iomodel->FetchMultipleData(&masscon_definitionstring_s,&nummasscons, "md.masscon.definitionstring");
451  iomodel->FetchMultipleData(&masscon_levelset_s,&masscon_levelset_M_s,&masscon_levelset_N_s,&nummasscons,"md.masscon.levelset");
452 
453  for(j=0;j<nummasscons;j++){
454 
455  /*Create a masscon object: */
456  output_definitions->AddObject(new Masscon(masscon_name_s[j],StringToEnumx(masscon_definitionstring_s[j]),masscon_levelset_s[j],masscon_levelset_M_s[j]));
457 
458  }
459 
460  /*Free ressources:*/
461  for(j=0;j<nummasscons;j++){
462  char* string=NULL;
463  IssmDouble* matrix = NULL;
464 
465  string = masscon_name_s[j]; xDelete<char>(string);
466  string = masscon_definitionstring_s[j]; xDelete<char>(string);
467  matrix = masscon_levelset_s[j]; xDelete<IssmDouble>(matrix);
468  }
469  xDelete<char*>(masscon_name_s);
470  xDelete<IssmDouble*>(masscon_levelset_s);
471  xDelete<int>(masscon_levelset_M_s);
472  xDelete<int>(masscon_levelset_N_s);
473  xDelete<char*>(masscon_definitionstring_s);
474 
475  /*}}}*/
476  }
477  else if (output_definition_enums[i]==MassconaxpbyEnum){
478  /*Deal with masscon combinations: {{{*/
479 
480  /*masscon variables: */
481  char** masscon_name_s = NULL;
482  char** masscon_definitionstring_s = NULL;
483  char** masscon_namex_s = NULL;
484  char** masscon_namey_s = NULL;
485  IssmDouble* masscon_alpha_s = NULL;
486  IssmDouble* masscon_beta_s = NULL;
487  int num;
488 
489  /*Fetch names and multiplicators, etc ... (see src/m/classes/masscon_axpby.m): */
490  iomodel->FetchMultipleData(&masscon_name_s,&num, "md.massconaxpby.name");
491  iomodel->FetchMultipleData(&masscon_definitionstring_s,&num,"md.massconaxpby.definitionstring");
492  iomodel->FetchMultipleData(&masscon_namex_s,&num, "md.massconaxpby.namex");
493  iomodel->FetchMultipleData(&masscon_namey_s,&num, "md.massconaxpby.namey");
494  iomodel->FetchMultipleData(&masscon_alpha_s,&num, "md.massconaxpby.alpha");
495  iomodel->FetchMultipleData(&masscon_beta_s,&num, "md.massconaxpby.beta");
496  for(j=0;j<num;j++){
497 
498  /*Create a masscon axpyb object: */
499  output_definitions->AddObject(new Massconaxpby(masscon_name_s[j],StringToEnumx(masscon_definitionstring_s[j]),masscon_namex_s[j],masscon_namey_s[j],masscon_alpha_s[j],masscon_beta_s[j]));
500 
501  }
502 
503  /*Free ressources:*/
504  for(j=0;j<num;j++){
505  char* string=NULL;
506  string = masscon_definitionstring_s[j]; xDelete<char>(string);
507  string = masscon_name_s[j]; xDelete<char>(string);
508  string = masscon_namex_s[j]; xDelete<char>(string);
509  string = masscon_namey_s[j]; xDelete<char>(string);
510  }
511  xDelete<char*>(masscon_definitionstring_s);
512  xDelete<char*>(masscon_name_s);
513  xDelete<char*>(masscon_namex_s);
514  xDelete<char*>(masscon_namey_s);
515  xDelete<IssmDouble>(masscon_alpha_s);
516  xDelete<IssmDouble>(masscon_beta_s);
517  /*}}}*/
518  }
519  else if (output_definition_enums[i]==RegionaloutputEnum){
520  /*Deal with regional output: {{{*/
521 
522  /*regional output variables: */
523  int numout;
524  char** reg_name_s = NULL;
525  char** reg_definitionstring_s = NULL;
526  char** reg_outputnamestring_s = NULL;
527  IssmDouble** reg_mask_s = NULL;
528  int* reg_mask_M_s = NULL;
529  int* reg_mask_N_s = NULL;
530 
531  /*Fetch name and mask, etc ... (see src/m/classes/regionaloutput.m): */
532  iomodel->FetchMultipleData(&reg_name_s,&numout, "md.regionaloutput.name");
533  iomodel->FetchMultipleData(&reg_definitionstring_s,&numout, "md.regionaloutput.definitionstring");
534  iomodel->FetchMultipleData(&reg_outputnamestring_s,&numout, "md.regionaloutput.outputnamestring");
535  iomodel->FetchMultipleData(&reg_mask_s,&reg_mask_M_s,&reg_mask_N_s,&numout, "md.regionaloutput.mask");
536  for(j=0;j<numout;j++){
537 
538  /*Create a regional output object: */
539  output_definitions->AddObject(new Regionaloutput(reg_name_s[j],StringToEnumx(reg_definitionstring_s[j]),reg_outputnamestring_s[j],reg_mask_s[j],reg_mask_M_s[j]));
540 
541  }
542 
543  /*Free ressources:*/
544  for(j=0;j<numout;j++){
545  char* string=NULL;
546  IssmDouble* matrix = NULL;
547 
548  string = reg_name_s[j]; xDelete<char>(string);
549  string = reg_definitionstring_s[j]; xDelete<char>(string);
550  string = reg_outputnamestring_s[j]; xDelete<char>(string);
551  matrix = reg_mask_s[j]; xDelete<IssmDouble>(matrix);
552  }
553  xDelete<char*>(reg_name_s);
554  xDelete<IssmDouble*>(reg_mask_s);
555  xDelete<int>(reg_mask_M_s);
556  xDelete<int>(reg_mask_N_s);
557  xDelete<char*>(reg_outputnamestring_s);
558  xDelete<char*>(reg_definitionstring_s);
559  /*}}}*/
560  }
561  else if (output_definition_enums[i]==NumberedcostfunctionEnum){
562  /*Deal with numbered cost function: {{{*/
563 
564  /*Intermediary*/
565  int numout,numout2;
566  char **ncf_name_s = NULL;
567  char **ncf_definitionstring_s = NULL;
568  char **cost_functions = NULL;
569  IssmDouble **cost_functions_weights = NULL;
570  int* cost_functions_weights_M = NULL;
571  int* cost_functions_weights_N = NULL;
572  int cost_function,domaintype;
573  int num_cost_functions;
574 
575  /*Process cost functions and convert from string to enums*/
576  iomodel->FindConstant(&num_cost_functions,"md.numberedcostfunction.num_cost_functions");
577  iomodel->FindConstant(&cost_functions,&num_cost_functions,"md.numberedcostfunction.cost_functions");
578  if(num_cost_functions<1) _error_("No cost functions found");
579  int* cost_function_enums=xNew<int>(num_cost_functions);
580  for(int i=0;i<num_cost_functions;++i){
581  cost_function_enums[i]=StringToEnumx(cost_functions[i]);
582  }
583 
584  iomodel->FetchMultipleData(&ncf_name_s,&numout,"md.numberedcostfunction.name");
585  iomodel->FetchMultipleData(&ncf_definitionstring_s,&numout2,"md.numberedcostfunction.definitionstring"); _assert_(numout2 == numout);
586  iomodel->FetchMultipleData(&cost_functions_weights,&cost_functions_weights_M,&cost_functions_weights_N,&numout2,"md.numberedcostfunction.cost_functions_coefficients"); _assert_(numout2 == numout);
587  if(numout!=1) _error_("not implemented yet, check code here");
588 
589  /*Fetch Observations */
590  iomodel->FindConstant(&domaintype,"md.mesh.domain_type");
591  for(int i=0;i<num_cost_functions;i++){
592  cost_function=cost_function_enums[i];
593  if( cost_function==ThicknessAbsMisfitEnum) iomodel->FetchDataToInput(inputs2,elements,"md.numberedcostfunction.thickness_obs",InversionThicknessObsEnum);
594  else if(cost_function==SurfaceAbsMisfitEnum) iomodel->FetchDataToInput(inputs2,elements,"md.numberedcostfunction.surface_obs",InversionSurfaceObsEnum);
595  else if(cost_function==SurfaceAbsVelMisfitEnum
596  || cost_function==SurfaceRelVelMisfitEnum
597  || cost_function==SurfaceLogVelMisfitEnum
598  || cost_function==SurfaceLogVxVyMisfitEnum
599  || cost_function==SurfaceAverageVelMisfitEnum){
600  iomodel->FetchDataToInput(inputs2,elements,"md.numberedcostfunction.vx_obs",InversionVxObsEnum);
601  if(domaintype!=Domain2DverticalEnum) iomodel->FetchDataToInput(inputs2,elements,"md.numberedcostfunction.vy_obs",InversionVyObsEnum);
602  }
603  }
604 
605  for(j=0;j<numout;j++){
606 
607  /*Now, for this particular misfit object, make sure we plug into the elements: the observation, and the weights.*/
608  for(int k=0;k<elements->Size();k++){
609  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(k));
610  element->DatasetInputCreate(cost_functions_weights[j],cost_functions_weights_M[j],cost_functions_weights_N[j],cost_function_enums,num_cost_functions,inputs2,iomodel,InversionCostFunctionsCoefficientsEnum);
611  }
612  output_definitions->AddObject(new Numberedcostfunction(ncf_name_s[j],StringToEnumx(ncf_definitionstring_s[j]),num_cost_functions,cost_function_enums));
613  }
614 
615  /*Free data: */
616  iomodel->DeleteData(2,"md.numberedcostfunction.name","md.numberedcostfunction.definitionstring");
617  xDelete<int>(cost_function_enums);
618  for(int i=0;i<num_cost_functions;i++) xDelete<char>(cost_functions[i]);
619  xDelete<char*>(cost_functions);
620 
621  /*Free ressources:*/
622  for(j=0;j<numout;j++){
623  xDelete<char>(ncf_name_s[j]);
624  xDelete<char>(ncf_definitionstring_s[j]);
625  xDelete<IssmDouble>(cost_functions_weights[j]);
626  }
627  xDelete<char*>(ncf_name_s);
628  xDelete<char*>(ncf_definitionstring_s);
629  xDelete<int>(cost_functions_weights_M);
630  xDelete<int>(cost_functions_weights_N);
631  xDelete<IssmDouble*>(cost_functions_weights);
632 
633  /*}}}*/
634  }
635  else if (output_definition_enums[i]==RadarEnum){
636  /*Deal with radar: {{{*/
637  int numout;
638  char **radar_name_s = NULL;
639  char **radar_definitionstring_s = NULL;
640  int **radar_ice_period_s = NULL;
641 
642  /*Fetch name and definition, etc ... (see src/m/classes/radar.m): */
643  iomodel->FetchMultipleData(&radar_definitionstring_s,&numout,"md.radar.definitionstring");
644  iomodel->FetchMultipleData(&radar_name_s,&numout,"md.radar.name");
645  if(numout>1) _error_("not suppored yet");
646  /*Fetch necessary inputs for calculation*/
647  //iomodel->FetchDataToInput(elements,"md.ice_period",RadarIcePeriodEnum);
648 
649  /*Add to output definitions*/
650  output_definitions->AddObject(new Radar(radar_name_s[0],StringToEnumx(radar_definitionstring_s[0])));
651  /*}}}*/
652  }
653  else _error_("output definition enum " << EnumToStringx(output_definition_enums[i]) << " not supported yet!");
654  }
655  }
656  parameters->AddObject(new DataSetParam(OutputdefinitionEnum,output_definitions));
657 
658  /*Free ressources:*/
659  delete output_definitions;
660  xDelete<int>(output_definition_enums);
661 
662 }
DataSet::Size
int Size()
Definition: DataSet.cpp:399
InversionThicknessObsEnum
@ InversionThicknessObsEnum
Definition: EnumDefinitions.h:631
Massfluxatgate
Definition: Massfluxatgate.h:18
SurfaceObservationEnum
@ SurfaceObservationEnum
Definition: EnumDefinitions.h:827
_assert_
#define _assert_(ignore)
Definition: exceptions.h:37
IssmDouble
double IssmDouble
Definition: types.h:37
MisfitEnum
@ MisfitEnum
Definition: EnumDefinitions.h:656
Cfdragcoeffabsgrad
Definition: Cfdragcoeffabsgrad.h:15
InversionVyObsEnum
@ InversionVyObsEnum
Definition: EnumDefinitions.h:634
RegionaloutputEnum
@ RegionaloutputEnum
Definition: EnumDefinitions.h:1237
DataSet::AddObject
int AddObject(Object *object)
Definition: DataSet.cpp:252
Nodalvalue
Definition: Nodalvalue.h:17
NodalvalueEnum
@ NodalvalueEnum
Definition: EnumDefinitions.h:1199
InversionCostFunctionsCoefficientsEnum
@ InversionCostFunctionsCoefficientsEnum
Definition: EnumDefinitions.h:629
CfdragcoeffabsgradEnum
@ CfdragcoeffabsgradEnum
Definition: EnumDefinitions.h:1005
SurfaceRelVelMisfitEnum
@ SurfaceRelVelMisfitEnum
Definition: EnumDefinitions.h:828
Radar
Definition: Radar.h:15
IoModel::FetchMultipleData
void FetchMultipleData(char ***pstringarray, int *pnumstrings, const char *data_name)
Definition: IoModel.cpp:1960
NumberedcostfunctionEnum
@ NumberedcostfunctionEnum
Definition: EnumDefinitions.h:1203
MassfluxatgateEnum
@ MassfluxatgateEnum
Definition: EnumDefinitions.h:1162
CfsurfacelogvelEnum
@ CfsurfacelogvelEnum
Definition: EnumDefinitions.h:1006
OutputdefinitionEnum
@ OutputdefinitionEnum
Definition: EnumDefinitions.h:285
Parameters::AddObject
void AddObject(Param *newparam)
Definition: Parameters.cpp:67
SurfaceLogVxVyMisfitEnum
@ SurfaceLogVxVyMisfitEnum
Definition: EnumDefinitions.h:826
SurfaceAbsVelMisfitEnum
@ SurfaceAbsVelMisfitEnum
Definition: EnumDefinitions.h:818
Element
Definition: Element.h:41
Numberedcostfunction
Definition: Numberedcostfunction.h:15
IoModel::numberofvertices
int numberofvertices
Definition: IoModel.h:99
MassconEnum
@ MassconEnum
Definition: EnumDefinitions.h:1160
Element::InputCreate
void InputCreate(IssmDouble *vector, Inputs2 *inputs2, IoModel *iomodel, int M, int N, int vector_type, int vector_enum, int code)
Definition: Element.cpp:1626
IoModel::DeleteData
void DeleteData(int num,...)
Definition: IoModel.cpp:500
IoModel::numberofelements
int numberofelements
Definition: IoModel.h:96
SurfaceAbsMisfitEnum
@ SurfaceAbsMisfitEnum
Definition: EnumDefinitions.h:817
VyObsEnum
@ VyObsEnum
Definition: EnumDefinitions.h:852
EnumToStringx
const char * EnumToStringx(int enum_in)
Definition: EnumToStringx.cpp:15
MassconaxpbyEnum
@ MassconaxpbyEnum
Definition: EnumDefinitions.h:1161
Cfsurfacesquare
Definition: Cfsurfacesquare.h:15
Regionaloutput
Definition: Regionaloutput.h:15
IoModel::FindConstant
void FindConstant(bool *pvalue, const char *constant_name)
Definition: IoModel.cpp:2362
ThicknessAbsMisfitEnum
@ ThicknessAbsMisfitEnum
Definition: EnumDefinitions.h:837
InversionSurfaceObsEnum
@ InversionSurfaceObsEnum
Definition: EnumDefinitions.h:630
IoModel::FetchData
void FetchData(bool *pboolean, const char *data_name)
Definition: IoModel.cpp:933
StringToEnumx
int StringToEnumx(const char *string_in, bool notfounderror=true)
Definition: StringToEnumx.cpp:14
DataSetParam
Definition: DataSetParam.h:20
Cfsurfacelogvel
Definition: Cfsurfacelogvel.h:15
RadarEnum
@ RadarEnum
Definition: EnumDefinitions.h:665
Element::DatasetInputCreate
virtual void DatasetInputCreate(IssmDouble *array, int M, int N, int *individual_enums, int num_inputs, Inputs2 *inputs2, IoModel *iomodel, int input_enum)
Definition: Element.h:218
Masscon
Definition: Masscon.h:18
_error_
#define _error_(StreamArgs)
Definition: exceptions.h:49
CfsurfacesquareEnum
@ CfsurfacesquareEnum
Definition: EnumDefinitions.h:1007
DataSet::GetObjectByOffset
Object * GetObjectByOffset(int offset)
Definition: DataSet.cpp:334
WeightsSurfaceObservationEnum
@ WeightsSurfaceObservationEnum
Definition: EnumDefinitions.h:864
SurfaceLogVelMisfitEnum
@ SurfaceLogVelMisfitEnum
Definition: EnumDefinitions.h:825
Massconaxpby
Definition: Massconaxpby.h:18
IoModel::FetchDataToInput
void FetchDataToInput(Inputs2 *inputs2, Elements *elements, const char *vector_name, int input_enum)
Definition: IoModel.cpp:1651
InversionVxObsEnum
@ InversionVxObsEnum
Definition: EnumDefinitions.h:633
SurfaceAverageVelMisfitEnum
@ SurfaceAverageVelMisfitEnum
Definition: EnumDefinitions.h:821
Element::DatasetInputAdd
void DatasetInputAdd(int enum_type, IssmDouble *vector, Inputs2 *inputs2, IoModel *iomodel, int M, int N, int vector_type, int vector_enum, int code, int input_enum)
Definition: Element.cpp:1831
DataSet
Declaration of DataSet class.
Definition: DataSet.h:14
Domain2DverticalEnum
@ Domain2DverticalEnum
Definition: EnumDefinitions.h:535
VxObsEnum
@ VxObsEnum
Definition: EnumDefinitions.h:848
Misfit
Definition: Misfit.h:15