Ice Sheet System Model  4.18
Code documentation
Public Member Functions | Data Fields | Private Attributes
IoModel Class Reference

#include <IoModel.h>

Public Member Functions

 ~IoModel ()
 
 IoModel ()
 
 IoModel (FILE *iomodel_handle, int solution_enum_in, bool trace, IssmPDouble *X)
 
void AddConstant (IoConstant *constant_in)
 
void AddConstantIndependent (IoConstant *constant_in)
 
void AddData (IoData *data_in)
 
void AddDataIndependent (IoData *data_in)
 
void FetchIndependentConstant (int *pXcount, IssmPDouble *X, const char *name)
 
void FetchIndependentData (int *pXcount, IssmPDouble *X, const char *name)
 
void FillIndependents (IssmDouble *xp)
 
void FindConstant (bool *pvalue, const char *constant_name)
 
void FindConstant (int *pvalue, const char *constant_name)
 
void FindConstant (IssmDouble *pvalue, const char *constant_name)
 
void FindConstant (char **pvalue, const char *constant_name)
 
void FindConstant (char ***pvalue, int *psize, const char *constant_name)
 
int NumIndependents ()
 
void CheckFile (void)
 
ParamCopyConstantObject (const char *constant_name, int param_enum)
 
IssmDoubleData (const char *data_name)
 
void DeclareIndependents (bool trace, IssmPDouble *X)
 
void DeleteData (int num,...)
 
void DeleteData (IssmDouble *vector, const char *data_name)
 
void DeleteData (char ***pstringarray, int numstrings, const char *data_name)
 
void FetchConstants (void)
 
void FetchData (bool *pboolean, const char *data_name)
 
void FetchData (int *pinteger, const char *data_name)
 
void FetchData (IssmDouble *pscalar, const char *data_name)
 
void FetchData (char **pstring, const char *data_name)
 
void FetchData (char ***pstrings, int *pnumstrings, const char *data_name)
 
void FetchData (int **pmatrix, int *pM, int *pN, const char *data_name)
 
void FetchData (IssmDouble **pscalarmatrix, int *pM, int *pN, const char *data_name)
 
void FetchData (IssmDouble ***pmatrixarray, int **pmdims, int **pndims, int *pnumrecords, const char *data_name)
 
void FetchData (Options *options, const char *data_name)
 
void FetchData (int num,...)
 
void FetchDataToInput (Inputs2 *inputs2, Elements *elements, const char *vector_name, int input_enum)
 
void FetchDataToInput (Inputs2 *inputs2, Elements *elements, const char *vector_name, int input_enum, IssmDouble default_value)
 
void FetchDataToDatasetInput (Inputs2 *inputs2, Elements *elements, const char *vector_name, int input_enum)
 
void FetchIndependent (const char *dependent_name)
 
void FetchMultipleData (char ***pstringarray, int *pnumstrings, const char *data_name)
 
void FetchMultipleData (IssmDouble ***pmatrixarray, int **pmdims, int **pndims, int *pnumrecords, const char *data_name)
 
void FetchMultipleData (int ***pmatrices, int **pmdims, int **pndims, int *pnumrecords, const char *data_name)
 
void FetchMultipleData (int **pvector, int *pnum_instances, const char *data_name)
 
void FetchMultipleData (IssmDouble **pvector, int *pnum_instances, const char *data_name)
 
fpos_t * SetFilePointersToData (int **pcodes, int **pvector_types, int *pnum_instances, const char *data_name)
 
FILE * SetFilePointerToData (int *pcode, int *pvector_type, const char *data_name)
 
void StartTrace (bool trace)
 

Data Fields

FILE * fid
 
int solution_enum
 
bool * my_elements
 
bool * my_faces
 
bool * my_vfaces
 
bool * my_edges
 
bool * my_vedges
 
bool * my_hedges
 
bool * my_vertices
 
int * my_vertices_lids
 
int * epart
 
int domaindim
 
int domaintype
 
int * elements
 
int * edges
 
int * verticaledges
 
int * horizontaledges
 
int * elementtoedgeconnectivity
 
int * elementtoverticaledgeconnectivity
 
int * elementtohorizontaledgeconnectivity
 
int * elementtofaceconnectivity
 
int * elementtoverticalfaceconnectivity
 
int * faces
 
int * verticalfaces
 
int facescols
 
int meshelementtype
 
int * numbernodetoelementconnectivity
 
int numberofedges
 
int numberofverticaledges
 
int numberofhorizontaledges
 
int numberofelements
 
int numberoffaces
 
int numberofverticalfaces
 
int numberofvertices
 
int * singlenodetoelementconnectivity
 

Private Attributes

std::vector< IoConstant * > constants
 
std::vector< IoData * > data
 

Detailed Description

Definition at line 48 of file IoModel.h.

Constructor & Destructor Documentation

◆ ~IoModel()

IoModel::~IoModel ( )

Definition at line 240 of file IoModel.cpp.

240  {/*{{{*/
241 
242  /*Delete constants*/
243  vector<IoConstant*>::iterator iter1;
244  for(iter1=constants.begin();iter1<constants.end();iter1++){
245  delete *iter1;
246  }
247  this->constants.clear();
248 
249  /*Delete data*/
250  vector<IoData*>::iterator iter2;
251  for(iter2=data.begin();iter2<data.end();iter2++){
252  #if defined(_ISSM_DEBUG_)
253  if(!(*iter2)->isindependent){
254  _printf0_("WARNING: IoData \"" << (*iter2)->name << "\" has not been freed (DeleteData has not been called)\n");
255  }
256  #endif
257  delete *iter2;
258  }
259  this->data.clear();
260 
261  xDelete<bool>(this->my_elements);
262  xDelete<bool>(this->my_faces);
263  xDelete<bool>(this->my_vfaces);
264  xDelete<bool>(this->my_edges);
265  xDelete<bool>(this->my_vedges);
266  xDelete<bool>(this->my_hedges);
267  xDelete<bool>(this->my_vertices);
268  xDelete<int>(this->my_vertices_lids);
269  xDelete<int>(this->epart);
270 
271  xDelete<int>(this->elements);
272  xDelete<int>(this->faces);
273  xDelete<int>(this->verticalfaces);
274  xDelete<int>(this->edges);
275  xDelete<int>(this->verticaledges);
276  xDelete<int>(this->horizontaledges);
277  xDelete<int>(this->elementtofaceconnectivity);
278  xDelete<int>(this->elementtoverticalfaceconnectivity);
279  xDelete<int>(this->elementtoedgeconnectivity);
280  xDelete<int>(this->elementtoverticaledgeconnectivity);
281  xDelete<int>(this->elementtohorizontaledgeconnectivity);
282  xDelete<int>(this->singlenodetoelementconnectivity);
283  xDelete<int>(this->numbernodetoelementconnectivity);
284 }

◆ IoModel() [1/2]

IoModel::IoModel ( )

Definition at line 131 of file IoModel.cpp.

131  {/*{{{*/
132 
133  this->fid=NULL;
134  this->solution_enum=-1;
135 
136  this->my_elements=NULL;
137  this->my_faces=NULL;
138  this->my_vfaces=NULL;
139  this->my_edges=NULL;
140  this->my_vedges=NULL;
141  this->my_hedges=NULL;
142  this->my_vertices=NULL;
143  this->my_vertices_lids=NULL;
144  this->epart=NULL;
145 
146  this->domaintype=-1;
147  this->domaindim=-1;
148  this->meshelementtype=-1;
149  this->numberofvertices=-1;
150  this->numberofelements=-1;
151  this->numberoffaces=-1;
152  this->numberofverticalfaces=-1;
153  this->numberofedges=-1;
154  this->numberofverticaledges=-1;
155  this->numberofhorizontaledges=-1;
156  this->facescols=-1;
157  this->elements=NULL;
158  this->faces=NULL;
159  this->verticalfaces=NULL;
160  this->edges=NULL;
161  this->verticaledges=NULL;
162  this->horizontaledges=NULL;
163  this->elementtofaceconnectivity = NULL;
165  this->elementtoedgeconnectivity = NULL;
168  this->singlenodetoelementconnectivity = NULL;
169  this->numbernodetoelementconnectivity = NULL;
170 }/*}}}*/

◆ IoModel() [2/2]

IoModel::IoModel ( FILE *  iomodel_handle,
int  solution_enum_in,
bool  trace,
IssmPDouble X 
)

Definition at line 171 of file IoModel.cpp.

171  {/*{{{*/
172 
173  bool autodiff=false;
174  bool iscontrol=false;
175 
176  /*First, keep track of the file handle: */
177  this->fid=iomodel_handle;
178 
179  /*Check that Enums are Synchronized*/
180  this->CheckFile();
181 
182  /*Keep track of solution*/
183  this->solution_enum = solution_enum_in;
184 
185  /*If we are running in AD mode, we need to start the trace and declare our independent variables now,
186  *and prevent them from being erased during successive calls to iomodel->FetchConstants, iomodel->FetchData and
187  iomodel->DeleteData:*/
188  this->StartTrace(trace);
189  this->DeclareIndependents(trace,X);
190 
191  /*Initialize and read constants:*/
192  this->FetchConstants(); /*this routine goes through the input file, and fetches bool, int, IssmDouble and string only, nothing memory intensive*/
193 
194  /*Is this an autodiff run?*/
195  this->FindConstant(&autodiff,"md.autodiff.isautodiff");
196  this->FindConstant(&iscontrol,"md.inversion.iscontrol");
197  if(trace){
198  autodiff=true;
199  }
200  else{
201  if(autodiff && !iscontrol)
202  autodiff=true;
203  else
204  autodiff=false;
205  }
206  this->AddConstant(new IoConstant(autodiff,"md.autodiff.isautodiff"));
207 
208  /*Initialize permanent data: */
209  this->my_elements = NULL;
210  this->my_faces = NULL;
211  this->my_vfaces = NULL;
212  this->my_edges = NULL;
213  this->my_vedges = NULL;
214  this->my_hedges = NULL;
215  this->my_vertices = NULL;
216  this->my_vertices_lids = NULL;
217  this->epart = NULL;
218 
219  FindConstant(&this->domaintype,"md.mesh.domain_type");
220  FindConstant(&this->meshelementtype,"md.mesh.elementtype");
221 
222  FetchData(&this->domaindim,"md.mesh.domain_dimension");
223  FetchData(&this->numberofvertices,"md.mesh.numberofvertices");
224  FetchData(&this->numberofelements,"md.mesh.numberofelements");
225  FetchData(&this->elements,NULL,NULL,"md.mesh.elements");
226  this->facescols = -1;
227  this->faces = NULL;
228  this->verticalfaces = NULL;
229  this->edges = NULL;
230  this->verticaledges = NULL;
231  this->horizontaledges = NULL;
232  this->elementtofaceconnectivity = NULL;
234  this->elementtoedgeconnectivity = NULL;
237  this->singlenodetoelementconnectivity = NULL;
238  this->numbernodetoelementconnectivity = NULL;
239 }/*}}}*/

Member Function Documentation

◆ AddConstant()

void IoModel::AddConstant ( IoConstant constant_in)

Definition at line 288 of file IoModel.cpp.

288  {/*{{{*/
289 
290  _assert_(in_constant);
291 
292  /*Go through dataset of constant and check whether it already exists */
293  vector<IoConstant*>::iterator iter;
294 
295  for(iter=constants.begin();iter<constants.end();iter++){
296  if(strcmp((*iter)->name,in_constant->name)==0){
297  delete in_constant;
298  return;
299  }
300  }
301 
302  this->constants.push_back(in_constant);
303 }

◆ AddConstantIndependent()

void IoModel::AddConstantIndependent ( IoConstant constant_in)

Definition at line 305 of file IoModel.cpp.

305  {/*{{{*/
306 
307  _assert_(in_constant);
308 
309  /*Set constant as independent*/
310  in_constant->isindependent = true;
311 
312  /*Add to constnats*/
313  this->AddConstant(in_constant);
314 }

◆ AddData()

void IoModel::AddData ( IoData data_in)

Definition at line 316 of file IoModel.cpp.

316  {/*{{{*/
317 
318  _assert_(in_data);
319 
320  /*Go through dataset of data and check whether it already exists */
321  vector<IoData*>::iterator iter;
322 
323  for(iter=data.begin();iter<data.end();iter++){
324  if(strcmp((*iter)->name,in_data->name)==0){
325  delete in_data;
326  return;
327  }
328  }
329 
330  this->data.push_back(in_data);
331 }

◆ AddDataIndependent()

void IoModel::AddDataIndependent ( IoData data_in)

Definition at line 333 of file IoModel.cpp.

333  {/*{{{*/
334 
335  _assert_(in_data);
336 
337  /*Set data as independent*/
338  in_data->isindependent = true;
339 
340  /*Add to constnats*/
341  this->AddData(in_data);
342 }

◆ FetchIndependentConstant()

void IoModel::FetchIndependentConstant ( int *  pXcount,
IssmPDouble X,
const char *  name 
)

Definition at line 1814 of file IoModel.cpp.

1814  {/*{{{*/
1815 
1816  /*recover my_rank:*/
1817  int my_rank=IssmComm::GetRank();
1818 
1819  /*recover Xcount if X is not NULL:*/
1820  int Xcount = 0;
1821  if(X) Xcount=*pXcount;
1822 
1823  #ifdef _HAVE_AD_ //cannot come here unless you are running AD mode, from DeclaredIndependents:
1824 
1825  /*output: */
1826  IssmPDouble pscalar;
1827  IssmDouble scalar; //same as pscalar, except it's an ADOLC independent variable
1828  int code;
1829 
1830  /*Set file pointer to beginning of the data: */
1831  fid=this->SetFilePointerToData(&code,NULL,constant_name);
1832  if(code!=3) _error_("expecting a IssmDouble for \"" << constant_name<<"\"");
1833 
1834  /*We have to read a scalar from disk. First read the dimensions of the scalar, then the scalar: */
1835  if(my_rank==0){
1836  if(fread(&pscalar,sizeof(IssmPDouble),1,fid)!=1)_error_("could not read scalar ");
1837 
1838  /*Now, before we even broadcast this to other nodes, declare the scalar as an independent variable!. If we
1839  *have been supplied an X vector, use it instead of what we just read: */
1840  #if defined(_HAVE_CODIPACK_)
1841  // FIXME codi here we just assign instead of using "operator <<="
1842  if(X){
1843  scalar=X[Xcount];
1844  } else {
1845  scalar=pscalar;
1846  }
1847  auto& tape_codi = IssmDouble::getGlobalTape();
1848  tape_codi.registerInput(scalar);
1849  codi_global.input_indices.push_back(scalar.getGradientData());
1850  #else
1851  if(X){
1852  scalar<<=X[Xcount];
1853  }
1854  else{
1855  scalar<<=pscalar;
1856  }
1857  #endif
1858  }
1859 
1861  this->AddConstantIndependent(new IoConstant(scalar,constant_name));
1862 
1863  /*increment offset into X vector, now that we have read 1 value:*/
1864  Xcount++; *pXcount=Xcount;
1865  #endif
1866 }

◆ FetchIndependentData()

void IoModel::FetchIndependentData ( int *  pXcount,
IssmPDouble X,
const char *  name 
)

Definition at line 1868 of file IoModel.cpp.

1868  {/*{{{*/
1869 
1870  /*recover my_rank:*/
1871  int my_rank=IssmComm::GetRank();
1872 
1873  /*recover Xcount if X is not NULL:*/
1874  int Xcount = 0;
1875  if(X) Xcount=*pXcount;
1876 
1877  #ifdef _HAVE_AD_ //cannot come here unless you are running AD mode, from DeclaredIndependents:
1878 
1879  /*Intermediaries*/
1880  int M,N;
1881  IssmPDouble* buffer=NULL; //a buffer to read the data from disk
1882  IssmDouble* matrix=NULL; //our independent variable
1883  int code,layout;
1884 
1885  /*Set file pointer to beginning of the data: */
1886  fid=this->SetFilePointerToData(&code,&layout,data_name);
1887  if((code!=5) && (code!=6) && (code!=7))_error_("expecting a IssmDouble, integer or boolean matrix for \"" << data_name<<"\"");
1888 
1889  /*We have to read a matrix from disk. First read the dimensions of the matrix, then the whole matrix: */
1890  /*numberofelements: */
1891  if(my_rank==0){
1892  if(fread(&M,sizeof(int),1,fid)!=1) _error_("could not read number of rows for matrix ");
1893  }
1895 
1896  if(my_rank==0){
1897  if(fread(&N,sizeof(int),1,fid)!=1) _error_("could not read number of columns for matrix ");
1898  }
1900 
1901  /*Now allocate matrix: */
1902  if(M*N){
1903  buffer=xNew<IssmPDouble>(M*N);
1904 // AD performance is sensitive to calls to ensurecontiguous.
1905 // Providing "t" will cause ensurecontiguous to be called.
1906 #ifdef _HAVE_AD_
1907  matrix=xNew<IssmDouble>(M*N,"t");
1908 #else
1909  matrix=xNew<IssmDouble>(M*N);
1910 #endif
1911 
1912  /*Read matrix on node 0, then broadcast: */
1913  if(my_rank==0){
1914  if(fread(buffer,M*N*sizeof(IssmPDouble),1,fid)!=1) _error_("could not read matrix ");
1915 
1916  /*Now, before we even broadcast this to other nodes, declare the whole matrix as a independent variable!
1917  If we have been supplied an X vector, use it instead of what we just read: */
1918  #if defined(_HAVE_CODIPACK_)
1919  // FIXME codi here we just assign instead of using "operator <<="
1920  auto& tape_codi = IssmDouble::getGlobalTape();
1921  if(X){
1922  for (int i=0;i<M*N;i++) {
1923  matrix[i]=X[Xcount+i];
1924  tape_codi.registerInput(matrix[i]);
1925  codi_global.input_indices.push_back(matrix[i].getGradientData());
1926  }
1927  }
1928  else{
1929  for (int i=0;i<M*N;i++) {
1930  matrix[i]=buffer[i];
1931  tape_codi.registerInput(matrix[i]);
1932  codi_global.input_indices.push_back(matrix[i].getGradientData());
1933  }
1934  }
1935  #else /*ADOLC*/
1936  if(X){
1937  for(int i=0;i<M*N;i++) matrix[i]<<=X[Xcount+i]; /*<<= ADOLC overloaded operator to declare independent*/
1938  }
1939  else{
1940  for(int i=0;i<M*N;i++) matrix[i]<<=buffer[i];
1941  }
1942  #endif
1943  }
1945 
1946  xDelete<IssmPDouble>(buffer);
1947  }
1948  else _error_("cannot declare the independent variable \"" << data_name << "\" if it's empty!");
1949 
1950  // FIXME codi is that at all relevant to CoDiPack or can we simply assume the same?
1951 
1952  /*Add to data as independent*/
1953  this->AddDataIndependent(new IoData(matrix,code,layout,M,N,data_name));
1954 
1955  /*increment offset into X vector, now that we have read M*N values:*/
1956  Xcount+=M*N; *pXcount=Xcount;
1957  #endif
1958 }

◆ FillIndependents()

void IoModel::FillIndependents ( IssmDouble xp)

Definition at line 2334 of file IoModel.cpp.

2334  {/*{{{*/
2335 
2336  _assert_(xp);
2337 
2338  /*Initialize local num ind*/
2339  int local_num_ind = 0;
2340 
2341  /*Process constants*/
2342  for(vector<IoConstant*>::iterator iter=constants.begin();iter<constants.end();iter++){
2343  if((*iter)->isindependent){
2344  (*iter)->constant->GetParameterValue(&xp[local_num_ind]);
2345  local_num_ind += 1;
2346  }
2347  }
2348 
2349  /*Process data*/
2350  for(vector<IoData*>::iterator iter=data.begin();iter<data.end();iter++){
2351  if((*iter)->isindependent){
2352  for(int i=0;i<(*iter)->M*(*iter)->N;i++){
2353  xp[local_num_ind+i] = (*iter)->data[i];
2354  }
2355  local_num_ind += (*iter)->M*(*iter)->N;
2356  }
2357  }
2358 
2359  _assert_(local_num_ind == this->NumIndependents());
2360 }

◆ FindConstant() [1/5]

void IoModel::FindConstant ( bool *  pvalue,
const char *  constant_name 
)

Definition at line 2362 of file IoModel.cpp.

2362  {/*{{{*/
2363 
2364  /*Intermediary*/
2365  vector<IoConstant*>::iterator iter;
2366 
2367  for(iter=constants.begin();iter<constants.end();iter++){
2368  IoConstant* ioconstant=*iter;
2369 
2370  if(strcmp(ioconstant->name,constant_name)==0){
2371  ioconstant->constant->GetParameterValue(pvalue);
2372  return;
2373  }
2374  }
2375 
2376  for(vector<IoConstant*>::iterator iter=constants.begin();iter<constants.end();iter++) (*iter)->constant->Echo();
2377  _error_("Could not find constant \""<<constant_name<<"\"");
2378 }

◆ FindConstant() [2/5]

void IoModel::FindConstant ( int *  pvalue,
const char *  constant_name 
)

Definition at line 2380 of file IoModel.cpp.

2380  {/*{{{*/
2381 
2382  /*Intermediary*/
2383  vector<IoConstant*>::iterator iter;
2384 
2385  for(iter=constants.begin();iter<constants.end();iter++){
2386  IoConstant* ioconstant=*iter;
2387 
2388  if(strcmp(ioconstant->name,constant_name)==0){
2389  ioconstant->constant->GetParameterValue(pvalue);
2390  return;
2391  }
2392  }
2393 
2394  _error_("Could not find constant \""<<constant_name <<"\"");
2395 }

◆ FindConstant() [3/5]

void IoModel::FindConstant ( IssmDouble pvalue,
const char *  constant_name 
)

Definition at line 2397 of file IoModel.cpp.

2397  {/*{{{*/
2398 
2399  /*Intermediary*/
2400  vector<IoConstant*>::iterator iter;
2401 
2402  for(iter=constants.begin();iter<constants.end();iter++){
2403  IoConstant* ioconstant=*iter;
2404 
2405  if(strcmp(ioconstant->name,constant_name)==0){
2406  ioconstant->constant->GetParameterValue(pvalue);
2407  return;
2408  }
2409  }
2410 
2411  _error_("Could not find constant \""<<constant_name <<"\"");
2412 }

◆ FindConstant() [4/5]

void IoModel::FindConstant ( char **  pvalue,
const char *  constant_name 
)

Definition at line 2414 of file IoModel.cpp.

2414  {/*{{{*/
2415 
2416  /*Intermediary*/
2417  vector<IoConstant*>::iterator iter;
2418 
2419  for(iter=constants.begin();iter<constants.end();iter++){
2420  IoConstant* ioconstant=*iter;
2421 
2422  if(strcmp(ioconstant->name,constant_name)==0){
2423  ioconstant->constant->GetParameterValue(pvalue);
2424  return;
2425  }
2426  }
2427 
2428  _error_("Could not find constant \""<<constant_name <<"\"");
2429 }

◆ FindConstant() [5/5]

void IoModel::FindConstant ( char ***  pvalue,
int *  psize,
const char *  constant_name 
)

Definition at line 2431 of file IoModel.cpp.

2431  {/*{{{*/
2432 
2433  /*Intermediary*/
2434  vector<IoConstant*>::iterator iter;
2435 
2436  for(iter=constants.begin();iter<constants.end();iter++){
2437  IoConstant* ioconstant=*iter;
2438 
2439  if(strcmp(ioconstant->name,constant_name)==0){
2440  ioconstant->constant->GetParameterValue(pvalue,psize);
2441  return;
2442  }
2443  }
2444 
2445  _error_("Could not find constant \""<<constant_name <<"\"");
2446 }

◆ NumIndependents()

int IoModel::NumIndependents ( void  )

Definition at line 2448 of file IoModel.cpp.

2448  {/*{{{*/
2449 
2450  /*Initialize output*/
2451  int num_independents = 0;
2452 
2453  /*Process constants*/
2454  for(vector<IoConstant*>::iterator iter=constants.begin();iter<constants.end();iter++){
2455  if((*iter)->isindependent){
2456  num_independents+= 1;
2457  }
2458  }
2459 
2460  /*Process data*/
2461  for(vector<IoData*>::iterator iter=data.begin();iter<data.end();iter++){
2462  if((*iter)->isindependent){
2463  num_independents+= (*iter)->M*(*iter)->N;
2464  }
2465  }
2466 
2467  /*return*/
2468  return num_independents;
2469 }

◆ CheckFile()

void IoModel::CheckFile ( void  )

Definition at line 344 of file IoModel.cpp.

344  {/*{{{*/
345 
346  bool found;
347  int record_enum,record_name_size;
348  long long record_length;
349  char *record_name = NULL;
350  const char *mddot = "md.";
351 
352  /*recover my_rank:*/
353  int my_rank=IssmComm::GetRank();
354 
355  /*Check that some fields have been allocated*/
356  _assert_(this->fid || my_rank);
357 
358  /*Go find in the binary file, the position of the data we want to fetch: */
359  if(my_rank==0){ //cpu 0
360 
361  /*First set FILE* position to the beginning of the file: */
362  fseek(this->fid,0,SEEK_SET);
363 
364  for(;;){
365  /*Read size of first string name: */
366  if(fread(&record_name_size,sizeof(int),1,fid)==0){
367  /*we have reached the end of the file. break: */
368  xDelete<char>(record_name);
369  break;
370  }
371  if(record_name_size<3 || record_name_size>80){
372  _error_("error while looking in binary file. Found a string of size "<<record_name_size);
373  }
374 
375  /*Allocate string of correct size: */
376  record_name=xNew<char>(record_name_size+1);
377  record_name[record_name_size]='\0';
378 
379  /*Read record_name: */
380  if(fread(record_name,record_name_size*sizeof(char),1,fid)==0){
381  /*we have reached the end of the file. break: */
382  found=false;
383  xDelete<char>(record_name);
384  break;
385  }
386  if(strncmp(record_name,mddot,3)!=0){
387  _error_("error while reading binary file: record does not start with \"md.\": "<<record_name);
388  }
389 
390  /*Have we found the last string?*/
391  if(strncmp(record_name,"md.EOF",6)==0){
392  found = true;
393  xDelete<char>(record_name);
394  break;
395  }
396 
397  /*Go to next Enum*/
398  if(fread(&record_length,sizeof(long long),1,fid)!=1) _error_("Could not read record_length");
399  fseek(fid,record_length,SEEK_CUR);
400  xDelete<char>(record_name);
401  }
402  if(!found){
403  _printf0_("\n");
404  _printf0_("=========================================================================\n");
405  _printf0_(" Marshalled file is corrupted \n");
406  _printf0_(" \n");
407  _printf0_(" * Last record found is : \n");
408  _printf0_(" the corresponding model field has probably been marshalled \n");
409  _printf0_(" incorrectly \n");
410  _printf0_(" \n");
411  _printf0_("=========================================================================\n");
412  _printf0_("\n");
413  _error_("Binary file corrupted (See error message above)");
414  }
415  }
416 }

◆ CopyConstantObject()

Param * IoModel::CopyConstantObject ( const char *  constant_name,
int  param_enum 
)

Definition at line 418 of file IoModel.cpp.

418  {/*{{{*/
419 
420  /*Intermediary*/
421  vector<IoConstant*>::iterator iter;
422 
423  for(iter=constants.begin();iter<constants.end();iter++){
424  IoConstant* ioconstant=*iter;
425 
426  if(strcmp(ioconstant->name,constant_name)==0){
427  Param* output = ioconstant->constant->copy();
428  output->SetEnum(param_enum);
429  return output;
430  }
431  }
432 
433  _error_("Constant \"" << constant_name << "\" not found in iomodel");
434  return NULL;
435 }

◆ Data()

IssmDouble * IoModel::Data ( const char *  data_name)

Definition at line 437 of file IoModel.cpp.

437  {/*{{{*/
438 
439  /*Intermediary*/
440  vector<IoData*>::iterator iter;
441 
442  for(iter=data.begin();iter<data.end();iter++){
443  IoData* iodata=*iter;
444  if(strcmp(iodata->name,data_name)==0) return iodata->data;
445  }
446 
447  return NULL;
448 }

◆ DeclareIndependents()

void IoModel::DeclareIndependents ( bool  trace,
IssmPDouble X 
)

Definition at line 450 of file IoModel.cpp.

450  {/*{{{*/
451 
452  bool autodiff,iscontrol;
453  int num_independent_objects,temp;
454  int Xcount=0;
455 
456  char** names = NULL;
457  int *types = NULL;
458 
459  /*Initialize array detecting whether data[i] is an independent AD mode variable: */
460  this->FetchData(&autodiff,"md.autodiff.isautodiff");
461  this->FetchData(&iscontrol,"md.inversion.iscontrol");
462 
463  if(trace || (autodiff && !iscontrol)){
464 
465  #ifdef _HAVE_AD_
466  // FIXME codi here we should be able to execute codi version as normal
467  this->FetchData(&num_independent_objects,"md.autodiff.num_independent_objects");
468  if(num_independent_objects){
469  this->FetchData(&names,&temp,"md.autodiff.independent_object_names");
470  _assert_(temp==num_independent_objects);
471  this->FetchData(&types,NULL,NULL,"md.autodiff.independent_object_types");
472 
473  /*create independent objects, and at the same time, fetch the corresponding independent variables,
474  *and declare them as such in ADOLC: */
475  for(int i=0;i<num_independent_objects;i++){
476 
477  if(types[i]==0){
478  /*Scalar*/
479  this->FetchIndependentConstant(&Xcount,X,names[i]);
480  }
481  else if(types[i]==1){
482  /* vector:*/
483  this->FetchIndependentData(&Xcount,X,names[i]);
484  }
485  else{
486  _error_("Independent cannot be of size " << types[i]);
487  }
488  }
489  for(int i=0;i<num_independent_objects;i++) xDelete<char>(names[i]);
490  xDelete<char*>(names);
491  xDelete<int>(types);
492  }
493  #else
494  /*if we asked for AD computations, we have a problem!: */
495  _error_("Cannot carry out AD mode computations without support of ADOLC or CoDiPack compiled in!");
496  #endif
497  }
498 }

◆ DeleteData() [1/3]

void IoModel::DeleteData ( int  num,
  ... 
)

Definition at line 500 of file IoModel.cpp.

500  {/*{{{*/
501 
502  /*Intermediaries*/
503  va_list ap;
504  char *data_name = NULL;
505  const char *mddot = "md.";
506  vector<IoData *>::iterator iter;
507 
508  /*Go through the entire list of data and delete the corresponding data from the iomodel-data dataset: */
509  va_start(ap,num);
510  for(int i=0;i<num;i++){
511  data_name=va_arg(ap,char*);
512 
513  if(strncmp(data_name,mddot,3)!=0) _error_("String provided does not start with \"md.\" ("<<data_name<<")");
514 
515  for(iter=data.begin();iter<data.end();iter++){
516  IoData* iodata=*iter;
517  if(strcmp(iodata->name,data_name)==0 && !iodata->isindependent){
518  delete *iter;
519  this->data.erase(iter);
520  break;
521  }
522  }
523  }
524  va_end(ap);
525 } /*}}}*/

◆ DeleteData() [2/3]

void IoModel::DeleteData ( IssmDouble vector,
const char *  data_name 
)

Definition at line 526 of file IoModel.cpp.

526  {/*{{{*/
527 
528  vector<IoData*>::iterator iter;
529 
530  /*do not do anything if pointer is NULL*/
531  if(!vector_in) return;
532 
533  /*do not delete if this is an independent variable*/
534  for(iter=data.begin();iter<data.end();iter++){
535  IoData* iodata=*iter;
536  if(strcmp(iodata->name,data_name)==0 && iodata->isindependent){
537  return;
538  }
539  }
540 
541  /*Go ahead and delete*/
542  xDelete<IssmDouble>(vector_in);
543 } /*}}}*/

◆ DeleteData() [3/3]

void IoModel::DeleteData ( char ***  pstringarray,
int  numstrings,
const char *  data_name 
)

Definition at line 544 of file IoModel.cpp.

544  {/*{{{*/
545 
546  char** stringarray=*pstringarray;
547 
548  if(numstrings){
549  for(int i=0;i<numstrings;i++){
550  char* string=stringarray[i];
551  xDelete<char>(string);
552  }
553  xDelete<char*>(stringarray);
554  }
555  *pstringarray=NULL;
556 } /*}}}*/

◆ FetchConstants()

void IoModel::FetchConstants ( void  )

Definition at line 557 of file IoModel.cpp.

557  {/*{{{*/
558 
559  /*record descriptions; */
560  const char* mddot = "md.";
561  char* record_name = NULL;
562  int record_name_size;
563  long long record_length;
564  int record_code; //1 to 7 number
565 
566  /*records: */
567  int booleanint = 0;
568  int integer = 0;
569  IssmPDouble pscalar = 0;
570  IssmDouble scalar = 0;
571  char *string = NULL;
572  char **strings = NULL;
573  int string_size,numstrings;
574 
575  /*recover my_rank:*/
576  int my_rank=IssmComm::GetRank();
577 
578  /*Check that some fields have been allocated*/
579  _assert_(this->fid || my_rank);
580 
581  /*Go find in the binary file, the position of the data we want to fetch: */
582  if(my_rank==0){ //cpu 0{{{
583 
584  /*First set FILE* position to the beginning of the file: */
585  fseek(this->fid,0,SEEK_SET);
586 
587  /*Now march through file looking for the correct data identifiers (bool,int,IssmDouble or string): */
588  for(;;){
589 
590  /*Read size of first string name: */
591  if(fread(&record_name_size,sizeof(int),1,fid)==0){
592  /*we have reached the end of the file. break: */
593  record_code=0; //0 means bailout
594  ISSM_MPI_Bcast(&record_code,1,ISSM_MPI_INT,0,IssmComm::GetComm()); /*tell others cpus we are bailing: */
595  break;
596  }
597  if(record_name_size<3 || record_name_size>80){
598  /*we are going to error out, still try and get informatoin for the user:*/
599  record_name=xNew<char>(record_name_size+1);
600  record_name[record_name_size]='\0';
601 
602  /*Read record_name: */
603  fread(record_name,record_name_size*sizeof(char),1,fid);
604 
605  _error_("error while looking in binary file. String " << record_name << " a string of size "<<record_name_size);
606  }
607 
608  /*Allocate string of correct size: */
609  record_name=xNew<char>(record_name_size+1);
610  record_name[record_name_size]='\0';
611 
612  /*Read record_name: */
613  if(fread(record_name,record_name_size*sizeof(char),1,fid)==0){
614  _error_("Could not read record name");
615  }
616 
617  if(strncmp(record_name,mddot,3)!=0){
618  _error_("error while reading binary file: record does not start with \"md.\": "<<record_name);
619  }
620 
621  /* Read the record length and the data type code: */
622  if(fread(&record_length,sizeof(long long),1,this->fid)!=1) _error_("Cound not read record_length");
623  if(fread(&record_code ,sizeof(int),1,this->fid)!=1) _error_("Cound not read record_code");
624 
625  /*Tell other cpus what we are doing: */
626  ISSM_MPI_Bcast(&record_code,1,ISSM_MPI_INT,0,IssmComm::GetComm()); /*tell other cpus what we are going to do: */
627 
628  /*Tell other cpus the name of the data, then branch according to the data type: */
629  ISSM_MPI_Bcast(&record_name_size,1,ISSM_MPI_INT,0,IssmComm::GetComm());
630  ISSM_MPI_Bcast(record_name,record_name_size,ISSM_MPI_CHAR,0,IssmComm::GetComm());
632 
633  switch(record_code){
634  case 1:
635  /*Read the boolean and broadcast it to other cpus:*/
636  if(fread(&booleanint,sizeof(int),1,this->fid)!=1) _error_("could not read boolean ");
637  ISSM_MPI_Bcast(&booleanint,1,ISSM_MPI_INT,0,IssmComm::GetComm());
638 
639  /*create BoolParam: */
640  this->AddConstant(new IoConstant((bool)booleanint,record_name)); //cast to boolean
641 
642  break;
643  case 2:
644  /*Read the integer and broadcast it to other cpus:*/
645  if(fread(&integer,sizeof(int),1,this->fid)!=1) _error_("could not read integer ");
646 
647  /*Convert codes to Enums if needed*/
648  if(strcmp(record_name,"md.smb.model")==0) integer = IoCodeToEnumSMB(integer);
649  if(strcmp(record_name,"md.basalforcings.model")==0) integer = IoCodeToEnumBasal(integer);
650  if(strcmp(record_name,"md.calving.law")==0) integer = IoCodeToEnumCalving(integer);
651  if(strcmp(record_name,"md.frontalforcings.parameterization")==0) integer = IoCodeToEnumFrontalforcings(integer);
652  if(strcmp(record_name,"md.hydrology.model")==0) integer = IoCodeToEnumHydrology(integer);
653  if(strcmp(record_name,"md.materials.type")==0) integer = IoCodeToEnumMaterials(integer);
654  if(strcmp(record_name,"md.materials.nature")==0) integer = IoCodeToEnumNature(integer);
655  if(strcmp(record_name,"md.timestepping.type")==0) integer = IoCodeToEnumTimestepping(integer);
656  if(strcmp(record_name,"md.amr.type")==0) integer = IoCodeToEnumAmr(integer);
657 
658  /*Broadcast to other cpus*/
660 
661  /*create IntParam: */
662  this->AddConstant(new IoConstant(integer,record_name));
663 
664  break;
665  case 3:
666  {
667  /*IssmDouble, check whether it is already there (from "declare independents")*/
668  bool exists = false;
669  vector<IoConstant*>::iterator iter;
670  for(iter=constants.begin();iter<constants.end();iter++){
671  IoConstant* ioconstant=*iter;
672  if(strcmp(ioconstant->name,record_name)==0){
673  exists = true;
674  break;
675  }
676  }
677  if(!exists){
678  if(fread(&pscalar,sizeof(IssmPDouble),1,this->fid)!=1) _error_("could not read scalar ");
680  scalar=pscalar;
681 
682  /*create DoubleParam: */
683  this->AddConstant(new IoConstant(scalar,record_name));
684  }
685  }
686  break;
687  case 4:
688  /*We have to read a string from disk. First read the dimensions of the string, then the string: */
689  if(fread(&string_size,sizeof(int),1,this->fid)!=1) _error_("could not read length of string ");
690  ISSM_MPI_Bcast(&string_size,1,ISSM_MPI_INT,0,IssmComm::GetComm());
691 
692  if(string_size){
693  string=xNew<char>(string_size+1);
694  string[string_size]='\0';
695 
696  /*Read string, then broadcast: */
697  if(fread(string,string_size*sizeof(char),1,this->fid)!=1)_error_(" could not read string ");
698  ISSM_MPI_Bcast(string,string_size,ISSM_MPI_CHAR,0,IssmComm::GetComm());
699  }
700  else{
701  string=xNew<char>(1);
702  string[0]='\0';
703  }
704  /*Convert strings to enums if needed*/
705  if(strcmp(record_name,"md.flowequation.fe_SSA")==0){
706  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
707  } else if(strcmp(record_name,"md.flowequation.fe_HO")==0){
708  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
709  } else if(strcmp(record_name,"md.flowequation.fe_FS")==0){
710  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
711  } else if(strcmp(record_name,"md.thermal.fe")==0){
712  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
713  } else if(strcmp(record_name,"md.levelset.fe")==0){
714  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
715  } else if(strcmp(record_name,"md.groundingline.migration")==0){
716  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
717  } else if(strcmp(record_name,"md.groundingline.friction_interpolation")==0){
718  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
719  } else if(strcmp(record_name,"md.groundingline.melt_interpolation")==0){
720  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
721  } else if(strcmp(record_name,"md.masstransport.hydrostatic_adjustment")==0){
722  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
723  } else if(strcmp(record_name,"md.materials.rheology_law")==0){
724  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
725  } else if(strcmp(record_name,"md.damage.elementinterp")==0){
726  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
727  } else if(strcmp(record_name,"md.mesh.domain_type")==0){
728  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
729  } else if(strcmp(record_name,"md.mesh.elementtype")==0){
730  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
731  } else {
732  /*Add string to parameters: */
733  this->AddConstant(new IoConstant(string,record_name));
734  }
735 
736  /*Free string*/
737  xDelete<char>(string);
738  break;
739  case 5:
740  case 6:
741  case 7:
742  case 8:
743  case 10:
744  /*We are not interested in this record, too memory intensive. Skip it: */
745  /*skip: */
746  fseek(fid,-sizeof(int),SEEK_CUR); //backtrak 1 integer
747  fseek(fid,record_length,SEEK_CUR);
748  break;
749  case 9:
750  /*String Array*/
751  if(fread(&numstrings,sizeof(int),1,fid)!=1) _error_("could not read length of string array");
752  ISSM_MPI_Bcast(&numstrings,1,ISSM_MPI_INT,0,IssmComm::GetComm());
753  /*Now allocate string array: */
754  if(numstrings){
755  strings=xNew<char*>(numstrings);
756  for(int i=0;i<numstrings;i++)strings[i]=NULL;
757 
758  /*Go through strings, and read: */
759  for(int i=0;i<numstrings;i++){
760 
761  if(fread(&string_size,sizeof(int),1,fid)!=1) _error_("could not read length of string ");
762  ISSM_MPI_Bcast(&string_size,1,ISSM_MPI_INT,0,IssmComm::GetComm());
763  if(string_size){
764  string=xNew<char>((string_size+1));
765  string[string_size]='\0';
766  if(fread(string,string_size*sizeof(char),1,fid)!=1)_error_(" could not read string ");
767  ISSM_MPI_Bcast(string,string_size,ISSM_MPI_CHAR,0,IssmComm::GetComm());
768  }
769  else{
770  string=xNew<char>(1);
771  string[0]='\0';
772  }
773  strings[i]=string;
774  }
775  }
776 
777  /*Add strings to parameters: */
778  this->AddConstant(new IoConstant(strings,numstrings,record_name));
779 
780  /*Free string*/
781  for(int i=0;i<numstrings;i++) xDelete<char>(strings[i]);
782  xDelete<char*>(strings);
783  break;
784  default:
785  _error_("unknown record type:" << record_code);
786  break;
787  }
788  xDelete<char>(record_name);
789  }
790  } //}}}
791  else{ //cpu ~0 {{{
792  for(;;){ //wait on cpu 0
793  ISSM_MPI_Bcast(&record_code,1,ISSM_MPI_INT,0,IssmComm::GetComm()); /*get from cpu 0 what we are going to do: */
794  if(record_code==0){
795  break; //we are done, break from the loop
796  }
797  else{
798  ISSM_MPI_Bcast(&record_name_size,1,ISSM_MPI_INT,0,IssmComm::GetComm());
799  _assert_(record_name_size);
800  record_name=xNew<char>((record_name_size+1)); record_name[record_name_size]='\0';
801  ISSM_MPI_Bcast(record_name,record_name_size,ISSM_MPI_CHAR,0,IssmComm::GetComm());
803  switch(record_code){
804  case 1:
805  /*boolean. get it from cpu 0 */
806  ISSM_MPI_Bcast(&booleanint,1,ISSM_MPI_INT,0,IssmComm::GetComm());
807 
808  /*create BoolParam: */
809  this->AddConstant(new IoConstant((bool)booleanint,record_name)); //cast to a boolean
810  break;
811 
812  case 2:
813  /*integer. get it from cpu 0 */
815 
816  /*create IntParam: */
817  this->AddConstant(new IoConstant(integer,record_name));
818  break;
819  case 3:
820  /*scalar. get it from cpu 0 */
821  {
822  /*IssmDouble, check whether it is already there (from "declare independents")*/
823  bool exists = false;
824  vector<IoConstant*>::iterator iter;
825  for(iter=constants.begin();iter<constants.end();iter++){
826  IoConstant* ioconstant=*iter;
827  if(strcmp(ioconstant->name,record_name)==0){
828  exists = true;
829  break;
830  }
831  }
832  if(!exists){
834  scalar=pscalar;
835  /*create DoubleParam: */
836  this->AddConstant(new IoConstant(scalar,record_name));
837  }
838  }
839  break;
840  case 4:
841  ISSM_MPI_Bcast(&string_size,1,ISSM_MPI_INT,0,IssmComm::GetComm());
842  if(string_size){
843  string=xNew<char>((string_size+1));
844  string[string_size]='\0';
845 
846  /*Read string from cpu 0: */
847  ISSM_MPI_Bcast(string,string_size,ISSM_MPI_CHAR,0,IssmComm::GetComm());
848  }
849  else{
850  string=xNew<char>(1);
851  string[0]='\0';
852  }
853 
854  if(strcmp(record_name,"md.flowequation.fe_SSA")==0){
855  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
856  } else if(strcmp(record_name,"md.flowequation.fe_HO")==0){
857  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
858  } else if(strcmp(record_name,"md.flowequation.fe_FS")==0){
859  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
860  } else if(strcmp(record_name,"md.thermal.fe")==0){
861  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
862  } else if(strcmp(record_name,"md.levelset.fe")==0){
863  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
864  } else if(strcmp(record_name,"md.groundingline.migration")==0){
865  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
866  } else if(strcmp(record_name,"md.groundingline.friction_interpolation")==0){
867  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
868  } else if(strcmp(record_name,"md.groundingline.melt_interpolation")==0){
869  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
870  } else if(strcmp(record_name,"md.masstransport.hydrostatic_adjustment")==0){
871  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
872  } else if(strcmp(record_name,"md.materials.rheology_law")==0){
873  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
874  } else if(strcmp(record_name,"md.damage.elementinterp")==0){
875  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
876  } else if(strcmp(record_name,"md.mesh.domain_type")==0){
877  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
878  } else if(strcmp(record_name,"md.mesh.elementtype")==0){
879  this->AddConstant(new IoConstant(StringToEnumx(string),record_name));
880  } else {
881  /*Add string to parameters: */
882  this->AddConstant(new IoConstant(string,record_name));
883  }
884 
885  /*Free string*/
886  xDelete<char>(string);
887  break;
888  case 5: break; //do nothing. not interested in this type of data, which is memory intensive.
889  case 6: break; //do nothing. not interested in this type of data, which is memory intensive.
890  case 7: break; //do nothing. not interested in this type of data, which is memory intensive.
891  case 8: break; //do nothing. not interested in this type of data, which is memory intensive.
892  case 10: break; //do nothing. not interested in this type of data, which is memory intensive.
893  case 9:
894  ISSM_MPI_Bcast(&numstrings,1,ISSM_MPI_INT,0,IssmComm::GetComm());
895  /*Now allocate string array: */
896  if(numstrings){
897  strings=xNew<char*>(numstrings);
898  for(int i=0;i<numstrings;i++)strings[i]=NULL;
899 
900  /*Go through strings, and read: */
901  for(int i=0;i<numstrings;i++){
902 
903  ISSM_MPI_Bcast(&string_size,1,ISSM_MPI_INT,0,IssmComm::GetComm());
904  if(string_size){
905  string=xNew<char>((string_size+1));
906  string[string_size]='\0';
907  ISSM_MPI_Bcast(string,string_size,ISSM_MPI_CHAR,0,IssmComm::GetComm());
908  }
909  else{
910  string=xNew<char>(1);
911  string[0]='\0';
912  }
913  strings[i]=string;
914  }
915  }
916 
917  /*Add strings to parameters: */
918  this->AddConstant(new IoConstant(strings,numstrings,record_name));
919 
920  /*Free string*/
921  for(int i=0;i<numstrings;i++) xDelete<char>(strings[i]);
922  xDelete<char*>(strings);
923  break;
924  default:
925  _error_("unknown record type:" << record_code);
926  break;
927  }
928  xDelete<char>(record_name);
929  }
930  }
931  } //}}}
932 }/*}}}*/

◆ FetchData() [1/10]

void IoModel::FetchData ( bool *  pboolean,
const char *  data_name 
)

Definition at line 933 of file IoModel.cpp.

933  {/*{{{*/
934 
935  /*output: */
936  int booleanint;
937  int code;
938 
939  /*recover my_rank:*/
940  int my_rank=IssmComm::GetRank();
941 
942  /*Set file pointer to beginning of the data: */
943  fid=this->SetFilePointerToData(&code,NULL,data_name);
944 
945  if(code!=1)_error_("expecting a boolean for \"" << data_name<<"\"");
946 
947  /*We have to read a boolean from disk. */
948  if(my_rank==0){
949  if(fread(&booleanint,sizeof(int),1,fid)!=1) _error_("could not read boolean ");
950  }
951  ISSM_MPI_Bcast(&booleanint,1,ISSM_MPI_INT,0,IssmComm::GetComm());
952 
953  /*cast to bool: */
954  /*Assign output pointers: */
955  *pboolean=(bool)booleanint;
956 
957 }

◆ FetchData() [2/10]

void IoModel::FetchData ( int *  pinteger,
const char *  data_name 
)

Definition at line 959 of file IoModel.cpp.

959  {/*{{{*/
960 
961  /*output: */
962  int integer;
963  int code;
964 
965  /*recover my_rank:*/
966  int my_rank=IssmComm::GetRank();
967  /*Set file pointer to beginning of the data: */
968  fid=this->SetFilePointerToData(&code,NULL,data_name);
969 
970  if(code!=2)_error_("expecting an integer for \"" << data_name<<"\"");
971 
972  /*We have to read a integer from disk. First read the dimensions of the integer, then the integer: */
973  if(my_rank==0){
974  if(fread(&integer,sizeof(int),1,fid)!=1) _error_("could not read integer ");
975  }
976 
978 
979  /*Assign output pointers: */
980  *pinteger=integer;
981 }

◆ FetchData() [3/10]

void IoModel::FetchData ( IssmDouble pscalar,
const char *  data_name 
)

Definition at line 983 of file IoModel.cpp.

983  {/*{{{*/
984 
985  /*output: */
986  IssmPDouble scalar;
987  int code;
988 
989  /*recover my_rank:*/
990  int my_rank=IssmComm::GetRank();
991 
992  /*Set file pointer to beginning of the data: */
993  fid=this->SetFilePointerToData(&code,NULL,data_name);
994 
995  if(code!=3)_error_("expecting a IssmDouble for \""<<data_name<<"\"");
996 
997  /*We have to read a scalar from disk. First read the dimensions of the scalar, then the scalar: */
998  if(my_rank==0){
999  if(fread(&scalar,sizeof(IssmPDouble),1,fid)!=1)_error_("could not read scalar ");
1000  }
1002 
1003  /*Assign output pointers: */
1004  *pscalar=scalar;
1005 
1006 }

◆ FetchData() [4/10]

void IoModel::FetchData ( char **  pstring,
const char *  data_name 
)

Definition at line 1008 of file IoModel.cpp.

1008  {/*{{{*/
1009 
1010  /*output: */
1011  char* string=NULL;
1012  int string_size;
1013  int code=0;
1014 
1015  /*recover my_rank:*/
1016  int my_rank=IssmComm::GetRank();
1017 
1018  /*Set file pointer to beginning of the data: */
1019  fid=this->SetFilePointerToData(&code,NULL,data_name);
1020 
1021  if(code!=4)_error_("expecting a string for \""<<data_name<<"\"");
1022 
1023  /*Now fetch: */
1024 
1025  /*We have to read a string from disk. First read the dimensions of the string, then the string: */
1026  if(my_rank==0){
1027  if(fread(&string_size,sizeof(int),1,fid)!=1) _error_("could not read length of string ");
1028  }
1029 
1030  ISSM_MPI_Bcast(&string_size,1,ISSM_MPI_INT,0,IssmComm::GetComm());
1031 
1032  /*Now allocate string: */
1033  if(string_size){
1034  string=xNew<char>((string_size+1));
1035  string[string_size]='\0';
1036 
1037  /*Read string on node 0, then broadcast: */
1038  if(my_rank==0){
1039  if(fread(string,string_size*sizeof(char),1,fid)!=1)_error_(" could not read string ");
1040  }
1041  ISSM_MPI_Bcast(string,string_size,ISSM_MPI_CHAR,0,IssmComm::GetComm());
1042  }
1043  else{
1044  string=xNew<char>(1);
1045  string[0]='\0';
1046  }
1047 
1048  /*Assign output pointers: */
1049  *pstring=string;
1050 }

◆ FetchData() [5/10]

void IoModel::FetchData ( char ***  pstrings,
int *  pnumstrings,
const char *  data_name 
)

Definition at line 1052 of file IoModel.cpp.

1052  {/*{{{*/
1053 
1054  /*output: */
1055  char** strings = NULL;
1056  char* string = NULL;
1057  int numstrings;
1058  int string_size;
1059  int code;
1060 
1061  /*recover my_rank:*/
1062  int my_rank=IssmComm::GetRank();
1063 
1064  /*Set file pointer to beginning of the data: */
1065  fid=this->SetFilePointerToData(&code,NULL,data_name);
1066 
1067  if(code!=9)_error_("expecting a string array for \""<<data_name<<"\"");
1068 
1069  /*Now fetch: */
1070 
1071  if(my_rank==0){
1072  if(fread(&numstrings,sizeof(int),1,fid)!=1) _error_("could not read length of string array");
1073  }
1074  ISSM_MPI_Bcast(&numstrings,1,ISSM_MPI_INT,0,IssmComm::GetComm());
1075 
1076  /*Now allocate string array: */
1077  if(numstrings){
1078  strings=xNew<char*>(numstrings);
1079  for(int i=0;i<numstrings;i++) strings[i]=NULL;
1080 
1081  /*Go through strings, and read: */
1082  for(int i=0;i<numstrings;i++){
1083 
1084  if(my_rank==0){
1085  if(fread(&string_size,sizeof(int),1,fid)!=1) _error_("could not read length of string ");
1086  }
1087  ISSM_MPI_Bcast(&string_size,1,ISSM_MPI_INT,0,IssmComm::GetComm());
1088  if(string_size){
1089  string=xNew<char>((string_size+1));
1090  string[string_size]='\0';
1091  if(my_rank==0){
1092  if(fread(string,string_size*sizeof(char),1,fid)!=1)_error_(" could not read string ");
1093  }
1094  ISSM_MPI_Bcast(string,string_size,ISSM_MPI_CHAR,0,IssmComm::GetComm());
1095  }
1096  else{
1097  string=xNew<char>(1);
1098  string[0]='\0';
1099  }
1100  strings[i]=string;
1101  }
1102  }
1103 
1104  /*Assign output pointers: */
1105  *pstrings = strings;
1106  if(pnumstrings) *pnumstrings = numstrings;
1107 }

◆ FetchData() [6/10]

void IoModel::FetchData ( int **  pmatrix,
int *  pM,
int *  pN,
const char *  data_name 
)

Definition at line 1109 of file IoModel.cpp.

1109  {/*{{{*/
1110  int i,j;
1111 
1112  /*output: */
1113  int M,N;
1114  IssmPDouble* matrix=NULL;
1115  int* integer_matrix=NULL;
1116  int code=0;
1117 
1118  /*recover my_rank:*/
1119  int my_rank=IssmComm::GetRank();
1120 
1121  /*Set file pointer to beginning of the data: */
1122  fid=this->SetFilePointerToData(&code,NULL,data_name);
1123 
1124  if(code!=5 && code!=6 && code!=7)_error_("expecting a IssmDouble, integer or boolean matrix for \""<<data_name<<"\""<<" (Code is "<<code<<")");
1125 
1126  /*Now fetch: */
1127 
1128  /*We have to read a matrix from disk. First read the dimensions of the matrix, then the whole matrix: */
1129  /*numberofelements: */
1130  if(my_rank==0){
1131  if(fread(&M,sizeof(int),1,fid)!=1) _error_("could not read number of rows for matrix ");
1132  }
1133 
1135 
1136  if(my_rank==0){
1137  if(fread(&N,sizeof(int),1,fid)!=1) _error_("could not read number of columns for matrix ");
1138  }
1140 
1141  /*Now allocate matrix: */
1142  if(M*N){
1143  matrix=xNew<IssmPDouble>(M*N);
1144 
1145  /*Read matrix on node 0, then broadcast: */
1146  if(my_rank==0){
1147  if(fread(matrix,M*N*sizeof(IssmPDouble),1,fid)!=1) _error_("could not read matrix ");
1148  }
1149 
1151  }
1152 
1153  /*Now cast to integer: */
1154  if(M*N){
1155  integer_matrix=xNew<int>(M*N);
1156  for (i=0;i<M;i++){
1157  for (j=0;j<N;j++){
1158  integer_matrix[i*N+j]=(int)matrix[i*N+j];
1159  }
1160  }
1161  }
1162  else{
1163  integer_matrix=NULL;
1164  }
1165  /*Free ressources:*/
1166  xDelete<IssmPDouble>(matrix);
1167 
1168  /*Assign output pointers: */
1169  *pmatrix=integer_matrix;
1170  if (pM)*pM=M;
1171  if (pN)*pN=N;
1172 
1173 }

◆ FetchData() [7/10]

void IoModel::FetchData ( IssmDouble **  pscalarmatrix,
int *  pM,
int *  pN,
const char *  data_name 
)

Definition at line 1175 of file IoModel.cpp.

1175  {/*{{{*/
1176 
1177  /*First, look if has already been loaded (might be an independent variable)*/
1178  vector<IoData*>::iterator iter;
1179  for(iter=data.begin();iter<data.end();iter++){
1180  IoData* iodata=*iter;
1181  if(strcmp(iodata->name,data_name)==0){
1182  *pmatrix=iodata->data;
1183  if(pM) *pM=iodata->M;
1184  if(pN) *pN=iodata->N;
1185  return;
1186  }
1187  }
1188 
1189  /*output: */
1190  int M,N;
1191  IssmPDouble *matrix = NULL;
1192  int code = 0;
1193 
1194  /*recover my_rank:*/
1195  int my_rank=IssmComm::GetRank();
1196 
1197  /*Set file pointer to beginning of the data: */
1198  fid=this->SetFilePointerToData(&code,NULL,data_name);
1199  if(code!=5 && code!=6 && code!=7 && code!=10)_error_("expecting a IssmDouble, integer or boolean matrix for \""<<data_name<<"\""<<" (Code is "<<code<<")");
1200 
1201  /*Now fetch: */
1202 
1203  /*We have to read a matrix from disk. First read the dimensions of the matrix, then the whole matrix: */
1204  /*numberofelements: */
1205  if(my_rank==0){
1206  if(fread(&M,sizeof(int),1,fid)!=1) _error_("could not read number of rows for matrix ");
1207  }
1209 
1210  if(my_rank==0){
1211  if(fread(&N,sizeof(int),1,fid)!=1) _error_("could not read number of columns for matrix ");
1212  }
1214 
1215  /*Now allocate matrix: */
1216  if(M*N){
1217  if(code==10){
1218  /*Special case for Compressed mat*/
1219  IssmPDouble offset,range;
1220  if(my_rank==0) if(fread(&offset,sizeof(IssmPDouble),1,fid)!=1) _error_("could not read offset");
1222 
1223  if(my_rank==0) if(fread(&range,sizeof(IssmPDouble),1,fid)!=1) _error_("could not read range");
1225 
1226  *pmatrix=xNew<IssmDouble>(M*N);
1227 
1228  /*Read matrix*/
1229  uint8_t* rawmatrix=xNew<uint8_t>((M-1)*N);
1230  if(my_rank==0) if(fread(rawmatrix,(M-1)*N*sizeof(char),1,fid)!=1) _error_("could not read matrix ");
1231  ISSM_MPI_Bcast(rawmatrix,(M-1)*N,ISSM_MPI_CHAR,0,IssmComm::GetComm());
1232 
1233  for(int i=0;i<(M-1)*N;++i) (*pmatrix)[i]=offset+range*reCast<IssmDouble>(rawmatrix[i])/255.;
1234  xDelete<uint8_t>(rawmatrix);
1235 
1236  /*read time now*/
1237  IssmPDouble* timematrix=xNew<IssmPDouble>(N);
1238  if(my_rank==0) if(fread(timematrix,N*sizeof(IssmPDouble),1,fid)!=1) _error_("could not read time in compressed matrix");
1240 
1241  for(int i=0;i<N;++i) (*pmatrix)[(M-1)*N+i]=timematrix[i];
1242  xDelete<IssmPDouble>(timematrix);
1243 
1244  }
1245  else{
1246  /*Read matrix on node 0, then broadcast: */
1247  matrix=xNew<IssmPDouble>(M*N);
1248  if(my_rank==0) if(fread(matrix,M*N*sizeof(IssmPDouble),1,fid)!=1) _error_("could not read matrix ");
1250 
1251  *pmatrix=xNew<IssmDouble>(M*N);
1252  for(int i=0;i<M*N;++i) (*pmatrix)[i]=matrix[i];
1253  xDelete<IssmPDouble>(matrix);
1254  }
1255  }
1256  else{
1257  *pmatrix=NULL;
1258  }
1259  /*Assign output pointers: */
1260  if(pM) *pM=M;
1261  if(pN) *pN=N;
1262 }

◆ FetchData() [8/10]

void IoModel::FetchData ( IssmDouble ***  pmatrixarray,
int **  pmdims,
int **  pndims,
int *  pnumrecords,
const char *  data_name 
)

Definition at line 1264 of file IoModel.cpp.

1264  {/*{{{*/
1265 
1266  int i;
1267  /*output: */
1268  IssmDouble** matrices=NULL;
1269  int* mdims=NULL;
1270  int* ndims=NULL;
1271  int numrecords=0;
1272 
1273  /*intermediary: */
1274  int M, N;
1275  IssmPDouble *matrix = NULL;
1276  int code;
1277 
1278  /*recover my_rank:*/
1279  int my_rank=IssmComm::GetRank();
1280 
1281  /*Set file pointer to beginning of the data: */
1282  fid=this->SetFilePointerToData(&code,NULL,data_name);
1283  if(code!=8)_error_("expecting a IssmDouble mat array for \""<<data_name<<"\"");
1284 
1285  /*Now fetch: */
1286  if(my_rank==0){
1287  if(fread(&numrecords,sizeof(int),1,fid)!=1) _error_("could not read number of records in matrix array ");
1288  }
1289  ISSM_MPI_Bcast(&numrecords,1,ISSM_MPI_INT,0,IssmComm::GetComm());
1290 
1291  if(numrecords){
1292 
1293  /*Allocate matrices :*/
1294  matrices=xNew<IssmDouble*>(numrecords);
1295  mdims=xNew<int>(numrecords);
1296  ndims=xNew<int>(numrecords);
1297 
1298  for(i=0;i<numrecords;i++){
1299  matrices[i]=NULL;
1300  mdims[i]=0;
1301  ndims[i]=0;
1302  }
1303 
1304  /*Loop through records and fetch matrix: */
1305  for(i=0;i<numrecords;i++){
1306 
1307  if(my_rank==0){
1308  if(fread(&M,sizeof(int),1,fid)!=1) _error_("could not read number of rows in " << i << "th matrix of matrix array");
1309  }
1311 
1312  if(my_rank==0){
1313  if(fread(&N,sizeof(int),1,fid)!=1) _error_("could not read number of columns in " << i << "th matrix of matrix array");
1314  }
1316 
1317  /*Now allocate matrix: */
1318  if(M*N){
1319  matrix=xNew<IssmPDouble>(M*N);
1320 
1321  /*Read matrix on node 0, then broadcast: */
1322  if(my_rank==0){
1323  if(fread(matrix,M*N*sizeof(IssmPDouble),1,fid)!=1) _error_("could not read matrix ");
1324  }
1325 
1327  matrices[i]=xNew<IssmDouble>(M*N);
1328  for (int j=0;j<M*N;++j) {matrices[i][j]=matrix[j];}
1329  xDelete<IssmPDouble>(matrix);
1330  }
1331  else
1332  matrices[i]=NULL;
1333  /*Assign: */
1334  mdims[i]=M;
1335  ndims[i]=N;
1336  }
1337  }
1338 
1339  /*Assign output pointers: */
1340  *pmatrices=matrices;
1341  *pmdims=mdims;
1342  *pndims=ndims;
1343  *pnumrecords=numrecords;
1344 }

◆ FetchData() [9/10]

void IoModel::FetchData ( Options options,
const char *  data_name 
)

Definition at line 1346 of file IoModel.cpp.

1346  {/*{{{*/
1347 
1348  /*record descriptions; */
1349  const char* mddot = "md.";
1350  char* record_name = NULL;
1351  int record_name_size;
1352  long long record_length;
1353  int record_code;
1354 
1355  /*records: */
1356  IssmDouble scalar = 0;
1357  char *string = NULL;
1358  int string_size;
1359 
1360  /*recover my_rank:*/
1361  int my_rank=IssmComm::GetRank();
1362 
1363  /*Go find in the binary file, the position of the data we want to fetch: */
1364  if(my_rank==0){
1365  fseek(fid,0,SEEK_SET);
1366  for(;;){
1367  /*Read size of first string name: */
1368  if(fread(&record_name_size,sizeof(int),1,fid)==0) _error_("could not read record_name");
1369  if(record_name_size<3 || record_name_size>80) _error_("error while looking in binary file. Found a string of size "<<record_name_size);
1370 
1371  /*Allocate string of correct size: */
1372  record_name=xNew<char>(record_name_size+1);
1373  record_name[record_name_size]='\0';
1374 
1375  /*Read record_name: */
1376  if(fread(record_name,record_name_size*sizeof(char),1,fid)==0)_error_("Could not find field "<<lastnonoption);
1377  if(strncmp(record_name,mddot,3)!=0) _error_("error while reading binary file: record does not start with \"md.\": "<<record_name);
1378 
1379  /*Is this the record sought for? : */
1380  if(strcmp(record_name,lastnonoption)==0){
1381  if(fread(&record_length,sizeof(long long),1,fid)!=1) _error_("Could not read record_length");
1382  fseek(fid,record_length,SEEK_CUR);
1383  xDelete<char>(record_name);
1384  break;
1385  }
1386  else{
1387  if(fread(&record_length,sizeof(long long),1,fid)!=1) _error_("Could not read record_length");
1388  fseek(fid,record_length,SEEK_CUR);
1389  xDelete<char>(record_name);
1390  }
1391  }
1392  }
1393 
1394  /*Go find in the binary file, the position of the data we want to fetch: */
1395  if(my_rank==0){ //cpu 0{{{
1396 
1397  /*Now march through file looking for the correct data identifiers (bool,int,IssmDouble or string): */
1398  for(;;){
1399 
1400  /*Read size of first string name: */
1401  if(fread(&record_name_size,sizeof(int),1,fid)==0){
1402  /*we have reached the end of the file. break: */
1403  record_code=0; //0 means bailout
1404  ISSM_MPI_Bcast(&record_code,1,ISSM_MPI_INT,0,IssmComm::GetComm()); /*tell others cpus we are bailing: */
1405  break;
1406  }
1407  if(record_name_size<3 || record_name_size>80){
1408  _error_("error while looking in binary file. Found a string of size "<<record_name_size);
1409  }
1410 
1411  /*Allocate string of correct size: */
1412  record_name=xNew<char>(record_name_size+1);
1413  record_name[record_name_size]='\0';
1414 
1415  /*Read record_name: */
1416  if(fread(record_name,record_name_size*sizeof(char),1,fid)==0){
1417  _error_("Could not read record name");
1418  }
1419  if(strncmp(record_name,mddot,3)!=0){
1420  _error_("error while reading binary file: record does not start with \"md.\": "<<record_name);
1421  }
1422  if(strcmp(record_name,"md.EOF")==0){
1423  xDelete<char>(record_name);
1424  record_code=0; //0 means bailout
1425  ISSM_MPI_Bcast(&record_code,1,ISSM_MPI_INT,0,IssmComm::GetComm()); /*tell others cpus we are bailing: */
1426  break;
1427  }
1428 
1429  /* Read the record length and the data type code: */
1430  if(fread(&record_length,sizeof(long long),1,this->fid)!=1) _error_("Cound not read record_length");
1431  if(fread(&record_code ,sizeof(int),1,this->fid)!=1) _error_("Cound not read record_code");
1432 
1433  /*Tell other cpus what we are doing: */
1434  ISSM_MPI_Bcast(&record_code,1,ISSM_MPI_INT,0,IssmComm::GetComm()); /*tell other cpus what we are going to do: */
1435 
1436  /*Tell other cpus the name of the data, then branch according to the data type: */
1437  ISSM_MPI_Bcast(&record_name_size,1,ISSM_MPI_INT,0,IssmComm::GetComm());
1438  ISSM_MPI_Bcast(record_name,record_name_size,ISSM_MPI_CHAR,0,IssmComm::GetComm());
1440 
1441  switch(record_code){
1442  case 3:
1443  {
1444  if(fread(&scalar,sizeof(IssmPDouble),1,this->fid)!=1) _error_("could not read scalar ");
1447  char* optionname=xNew<char>(strlen(record_name)-3+1);
1448  xMemCpy(optionname,&record_name[3],strlen(record_name)-3+1);
1449  option->value = scalar;
1450  option->name = optionname;
1451  option->size[0] = 1;
1452  option->size[1] = 1;
1453  options->AddOption(option);
1454  }
1455  break;
1456  case 4:
1457  {
1458  /*We have to read a string from disk. First read the dimensions of the string, then the string: */
1459  if(fread(&string_size,sizeof(int),1,this->fid)!=1) _error_("could not read length of string ");
1460  ISSM_MPI_Bcast(&string_size,1,ISSM_MPI_INT,0,IssmComm::GetComm());
1461 
1462  if(string_size){
1463  string=xNew<char>(string_size+1);
1464  string[string_size]='\0';
1465 
1466  /*Read string, then broadcast: */
1467  if(fread(string,string_size*sizeof(char),1,this->fid)!=1)_error_(" could not read string ");
1468  ISSM_MPI_Bcast(string,string_size,ISSM_MPI_CHAR,0,IssmComm::GetComm());
1469  }
1470  else{
1471  string=xNew<char>(1);
1472  string[0]='\0';
1473  }
1474 
1475  /*Add string to parameters: */
1477  char* optionname=xNew<char>(strlen(record_name)-3+1);
1478  xMemCpy(optionname,&record_name[3],strlen(record_name)-3+1);
1479  option->value = string;
1480  option->name = optionname;
1481  option->size[0] = 1;
1482  option->size[1] = 1;
1483  options->AddOption(option);
1484 
1485  }
1486  break;
1487  default:
1488  _error_("record type not supported:" << record_code);
1489  break;
1490  }
1491  xDelete<char>(record_name);
1492  }
1493  } //}}}
1494  else{ //cpu ~0 {{{
1495  for(;;){ //wait on cpu 0
1496  ISSM_MPI_Bcast(&record_code,1,ISSM_MPI_INT,0,IssmComm::GetComm()); /*get from cpu 0 what we are going to do: */
1497  if(record_code==0){
1498  break; //we are done, break from the loop
1499  }
1500  else{
1501  ISSM_MPI_Bcast(&record_name_size,1,ISSM_MPI_INT,0,IssmComm::GetComm());
1502  _assert_(record_name_size);
1503  record_name=xNew<char>((record_name_size+1)); record_name[record_name_size]='\0';
1504  ISSM_MPI_Bcast(record_name,record_name_size,ISSM_MPI_CHAR,0,IssmComm::GetComm());
1506  switch(record_code){
1507  case 3:
1508  {
1509  if(fread(&scalar,sizeof(IssmPDouble),1,this->fid)!=1) _error_("could not read scalar ");
1511  char* optionname=xNew<char>(strlen(record_name)-3+1);
1512  xMemCpy(optionname,&record_name[3],strlen(record_name)-3+1);
1514  option->value = scalar;
1515  option->name = optionname;
1516  option->size[0] = 1;
1517  option->size[1] = 1;
1518  options->AddOption(option);
1519  }
1520  break;
1521  case 4:
1522  {
1523  /*We have to read a string from disk. First read the dimensions of the string, then the string: */
1524  if(fread(&string_size,sizeof(int),1,this->fid)!=1) _error_("could not read length of string ");
1525  ISSM_MPI_Bcast(&string_size,1,ISSM_MPI_INT,0,IssmComm::GetComm());
1526 
1527  if(string_size){
1528  string=xNew<char>(string_size+1);
1529  string[string_size]='\0';
1530 
1531  /*Read string, then broadcast: */
1532  if(fread(string,string_size*sizeof(char),1,this->fid)!=1)_error_(" could not read string ");
1533  ISSM_MPI_Bcast(string,string_size,ISSM_MPI_CHAR,0,IssmComm::GetComm());
1534  }
1535  else{
1536  string=xNew<char>(1);
1537  string[0]='\0';
1538  }
1539 
1540  /*Add string to parameters: */
1541  char* optionname=xNew<char>(strlen(record_name)-3+1);
1542  xMemCpy(optionname,&record_name[3],strlen(record_name)-3+1);
1544  option->value = string;
1545  option->name = optionname;
1546  option->size[0] = 1;
1547  option->size[1] = 1;
1548  options->AddOption(option);
1549  }
1550  break;
1551  default:
1552  _error_("record type not supported:" << record_code);
1553  break;
1554  }
1555 
1556  }
1557  }
1558  } //}}}
1559 
1560 }

◆ FetchData() [10/10]

void IoModel::FetchData ( int  num,
  ... 
)

Definition at line 1562 of file IoModel.cpp.

1562  {/*{{{*/
1563 
1564  va_list ap;
1565  int code,layout;
1566  IssmDouble *matrix = NULL;
1567  char* data_name;
1568  int M,N;
1569  bool exists;
1570  const char *mddot = "md.";
1571  vector<IoData*>::iterator iter;
1572 
1573  /*Go through the entire list of names and fetch the corresponding data. Add it to the iomodel->data dataset. Everything
1574  *we fetch is a IssmDouble* : */
1575 
1576  va_start(ap,num);
1577  for(int i=0; i<num; i++){
1578 
1579  data_name=va_arg(ap,char*);
1580  if(strncmp(data_name,mddot,3)!=0) _error_("String provided does not start with \"md.\" ("<<data_name<<")");
1581 
1582  exists = false;
1583 
1584  for(iter=data.begin();iter<data.end();iter++){
1585  IoData* iodata=*iter;
1586  if(strcmp(iodata->name,data_name)==0){
1587  /*Already there, no need to fetch it*/
1588  _assert_(iodata->isindependent);
1589  exists = true;
1590  break;
1591  }
1592  }
1593 
1594  if(exists){
1595  /*this data has already been checked out! Continue: */
1596  continue;
1597  }
1598  else{
1599  /*Add to this->data: */
1600  this->SetFilePointerToData(&code,&layout,data_name);
1601  this->FetchData(&matrix,&M,&N,data_name);
1602  this->AddData(new IoData(matrix,code,layout,M,N,data_name));
1603  }
1604  }
1605  va_end(ap);
1606 }

◆ FetchDataToInput() [1/2]

void IoModel::FetchDataToInput ( Inputs2 inputs2,
Elements elements,
const char *  vector_name,
int  input_enum 
)

Definition at line 1651 of file IoModel.cpp.

1651  {/*{{{*/
1652 
1653  /*First, look whether it is not already loaded in this->data*/
1654  vector<IoData*>::iterator iter;
1655  for(iter=data.begin();iter<data.end();iter++){
1656  IoData* iodata=*iter;
1657  if(strcmp(iodata->name,vector_name)==0){
1658  for(int i=0;i<elements->Size();i++){
1659  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(i));
1660  element->InputCreate(iodata->data,inputs2,this,iodata->M,iodata->N,iodata->layout,input_enum,iodata->code);//we need i to index into elements.
1661  }
1662  return;
1663  }
1664  }
1665 
1666  /*intermediary: */
1667  int i;
1668  int code,vector_layout;
1669 
1670  /*variables being fetched: */
1671  bool boolean;
1672  int integer;
1673  IssmDouble scalar;
1674  IssmDouble *doublearray = NULL;
1675  int M,N;
1676 
1677  /*First of, find the record for the name, and get code of data type: */
1678  this->SetFilePointerToData(&code, &vector_layout,vector_name);
1679 
1680  switch(code){
1681  case 1: //boolean constant
1682  this->FetchData(&boolean,vector_name);
1683  for(i=0;i<elements->Size();i++){
1684  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(i));
1685  element->SetBoolInput(inputs2,input_enum,boolean);
1686  }
1687  break;
1688  case 2: //integer constant
1689  this->FetchData(&integer,vector_name);
1690  for(i=0;i<elements->Size();i++){
1691  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(i));
1692  element->SetIntInput(inputs2,input_enum,integer);
1693  }
1694  break;
1695  case 3: //IssmDouble constant
1696  this->FetchData(&scalar,vector_name);
1697  for(i=0;i<elements->Size();i++){
1698  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(i));
1699  element->SetElementInput(inputs2,input_enum,scalar);
1700  }
1701  break;
1702  case 5: //boolean vector
1703  this->FetchData(&doublearray,&M,&N,vector_name); //we still have a doublearray, because it might include times in transient mode
1704  if(!doublearray) _error_("\""<<vector_name<<"\" not found in binary file");
1705  for(i=0;i<elements->Size();i++){
1706  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(i));
1707  element->InputCreate(doublearray,inputs2,this,M,N,vector_layout,input_enum,code);//we need i to index into elements.
1708  }
1709  break;
1710  case 6: //int vector
1711  this->FetchData(&doublearray,&M,&N,vector_name); //we still have a doublearray, because it might include times in transient mode
1712  if(!doublearray) _error_("\""<<vector_name<<"\" not found in binary file");
1713  for(i=0;i<elements->Size();i++){
1714  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(i));
1715  element->InputCreate(doublearray,inputs2,this,M,N,vector_layout,input_enum,code);//we need i to index into elements.
1716  }
1717  break;
1718  case 7: //IssmDouble vector
1719  case 10:
1720  this->FetchData(&doublearray,&M,&N,vector_name);
1721  if(!doublearray) _error_("\""<<vector_name<<"\" not found in binary file");
1722  for(i=0;i<elements->Size();i++){
1723  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(i));
1724  element->InputCreate(doublearray,inputs2,this,M,N,vector_layout,input_enum,code);//we need i to index into elements.
1725  }
1726  break;
1727  default:
1728  _error_("data code " << code << " not supported yet (detected while processing \""<<vector_name<<"\")");
1729  break;
1730  }
1731  /*Free ressources*/
1732  xDelete<IssmDouble>(doublearray);
1733 }

◆ FetchDataToInput() [2/2]

void IoModel::FetchDataToInput ( Inputs2 inputs2,
Elements elements,
const char *  vector_name,
int  input_enum,
IssmDouble  default_value 
)

Definition at line 1608 of file IoModel.cpp.

1608  {/*{{{*/
1609 
1610  /*First, look whether it is not already loaded in this->data*/
1611  vector<IoData*>::iterator iter;
1612  for(iter=data.begin();iter<data.end();iter++){
1613  IoData* iodata=*iter;
1614  if(strcmp(iodata->name,vector_name)==0){
1615  _assert_(iodata->code==7);
1616  for(int i=0;i<elements->Size();i++){
1617  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(i));
1618  element->InputCreate(iodata->data,inputs2,this,iodata->M,iodata->N,iodata->layout,input_enum,iodata->code);//we need i to index into elements.
1619  }
1620  return;
1621  }
1622  }
1623 
1624  /*intermediary: */
1625  int code,vector_layout;
1626  IssmDouble *doublearray = NULL;
1627  int M,N;
1628 
1629  /*First of, find the record for the name, and get code of data type: */
1630  this->SetFilePointerToData(&code, &vector_layout,vector_name);
1631 
1632  /*Defaulting only supported for double arrays*/
1633  if(code!=7 && code!=10) _error_(vector_name<<" is not a double array");
1634 
1635  this->FetchData(&doublearray,&M,&N,vector_name);
1636 
1637  for(int i=0;i<elements->Size();i++){
1638  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(i));
1639  if(!doublearray){
1640  element->SetElementInput(inputs2,input_enum,default_value);
1641  }
1642  else{
1643  element->InputCreate(doublearray,inputs2,this,M,N,vector_layout,input_enum,code);//we need i to index into elements.
1644  }
1645  }
1646 
1647  /*Free ressources*/
1648  xDelete<IssmDouble>(doublearray);
1649 }

◆ FetchDataToDatasetInput()

void IoModel::FetchDataToDatasetInput ( Inputs2 inputs2,
Elements elements,
const char *  vector_name,
int  input_enum 
)

Definition at line 1735 of file IoModel.cpp.

1735  {/*{{{*/
1736 
1737  /*First, look whether it is not already loaded in this->data*/
1738  vector<IoData*>::iterator iter;
1739  for(iter=data.begin();iter<data.end();iter++){
1740  IoData* iodata=*iter;
1741  if(strcmp(iodata->name,vector_name)==0){
1742  for(int i=0;i<elements->Size();i++){
1743  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(i));
1744  _error_("to be implemented...");
1745  //element->InputCreate(iodata->data,inputs2,this,iodata->M,iodata->N,iodata->layout,input_enum,iodata->code);//we need i to index into elements.
1746  }
1747  return;
1748  }
1749  }
1750 
1751  /*intermediary: */
1752  int code,vector_layout;
1753  IssmDouble *doublearray = NULL;
1754  int M,N;
1755 
1756  /*First of, find the record for the name, and get code of data type: */
1757  this->SetFilePointerToData(&code,&vector_layout,vector_name);
1758 
1759  switch(code){
1760  case 1: //boolean constant
1761  _error_("not implemented yet");
1762  break;
1763  case 2: //integer constant
1764  _error_("not implemented yet");
1765  break;
1766  case 3: //IssmDouble constant
1767  _error_("not implemented yet");
1768  break;
1769  case 5: //boolean vector
1770  _error_("not implemented yet");
1771  break;
1772  case 6: //int vector
1773  _error_("not implemented yet");
1774  break;
1775  case 7: //IssmDouble vector
1776  {
1777  this->FetchData(&doublearray,&M,&N,vector_name);
1778  if(!doublearray) _error_("\""<<vector_name<<"\" not found in binary file");
1779 
1780  int* ids = xNew<int>(N);
1781  for(int i=0;i<N;i++) ids[i] = i;
1782 
1783  for(int i=0;i<elements->Size();i++){
1784  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(i));
1785  element->DatasetInputCreate(doublearray,M,N,ids,N,inputs2,this,input_enum);
1786  }
1787  xDelete<int>(ids);
1788  }
1789  break;
1790  case 10: //Compressed matrix
1791  {
1792  this->FetchData(&doublearray,&M,&N,vector_name);
1793  if(!doublearray) _error_("\""<<vector_name<<"\" not found in binary file");
1794 
1795  int* ids = xNew<int>(N);
1796  for(int i=0;i<N;i++) ids[i] = i;
1797 
1798  for(int i=0;i<elements->Size();i++){
1799  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(i));
1800  element->DatasetInputCreate(doublearray,M,N,ids,N,inputs2,this,input_enum);
1801  }
1802  xDelete<int>(ids);
1803  }
1804  break;
1805 
1806  default:
1807  _error_("data code " << code << " not supported yet (detected while processing \""<<vector_name<<"\")");
1808  break;
1809  }
1810  /*Free ressources*/
1811  xDelete<IssmDouble>(doublearray);
1812 }

◆ FetchIndependent()

void IoModel::FetchIndependent ( const char *  dependent_name)

◆ FetchMultipleData() [1/5]

void IoModel::FetchMultipleData ( char ***  pstringarray,
int *  pnumstrings,
const char *  data_name 
)

Definition at line 1960 of file IoModel.cpp.

1960  {/*{{{*/
1961 
1962  int num_instances;
1963 
1964  /*output: */
1965  int numstrings = 0;
1966  char **strings = NULL;
1967 
1968  /*intermediary: */
1969  char *string = NULL;
1970  int string_size;
1971  int *codes = NULL;
1972  int *code = NULL;
1973  fpos_t *file_positions = NULL;
1974 
1975  /*recover my_rank:*/
1976  int my_rank=IssmComm::GetRank();
1977 
1978  /*Get file pointers to beginning of the data (multiple instances of it): */
1979  file_positions=this->SetFilePointersToData(&codes,NULL,&num_instances,data_name);
1980 
1981  if(num_instances){
1982  strings=xNew<char*>(num_instances);
1983 
1984  for(int i=0;i<num_instances;i++){
1985 
1986  if(my_rank==0){
1987  /*check we are indeed finding a string, not something else: */
1988  if(codes[i]!=4)_error_("expecting a string for \""<<data_name<<"\" but code is "<<codes[i]<<" not 4");
1989 
1990  /*We have to read a string from disk. First read the dimensions of the string, then the string: */
1991  fsetpos(fid,file_positions+i);
1992  if(fread(&string_size,sizeof(int),1,fid)!=1) _error_("could not read length of string ");
1993  }
1994 
1995  ISSM_MPI_Bcast(&string_size,1,ISSM_MPI_INT,0,IssmComm::GetComm());
1996 
1997  /*Now allocate string: */
1998  if(string_size){
1999  string=xNew<char>((string_size+1));
2000  string[string_size]='\0';
2001 
2002  /*Read string on node 0, then broadcast: */
2003  if(my_rank==0){
2004  if(fread(string,string_size*sizeof(char),1,fid)!=1)_error_(" could not read string ");
2005  }
2006  ISSM_MPI_Bcast(string,string_size,ISSM_MPI_CHAR,0,IssmComm::GetComm());
2007  }
2008  else{
2009  string=xNew<char>(1);
2010  string[0]='\0';
2011  }
2012  strings[i]=string;
2013  }
2014  }
2015  /*Free ressources:*/
2016  xDelete<int>(codes);
2017  xDelete<fpos_t>(file_positions);
2018 
2019  /*Assign output pointers: */
2020  *pstrings=strings;
2021  *pnumstrings=num_instances;
2022 }

◆ FetchMultipleData() [2/5]

void IoModel::FetchMultipleData ( IssmDouble ***  pmatrixarray,
int **  pmdims,
int **  pndims,
int *  pnumrecords,
const char *  data_name 
)

Definition at line 2130 of file IoModel.cpp.

2130  {/*{{{*/
2131 
2132  int num_instances;
2133  fpos_t* file_positions=NULL;
2134 
2135  /*output: */
2136  IssmDouble **matrices = NULL;
2137  int *mdims = NULL;
2138  int *ndims = NULL;
2139 
2140  /*intermediary: */
2141  int M, N;
2142  IssmPDouble *pmatrix = NULL;
2143  IssmDouble *matrix = NULL;
2144  int *codes = NULL;
2145  int code;
2146 
2147  /*recover my_rank:*/
2148  int my_rank=IssmComm::GetRank();
2149 
2150  /*Get file pointers to beginning of the data (multiple instances of it): */
2151  file_positions=this->SetFilePointersToData(&codes,NULL,&num_instances,data_name);
2152 
2153  if(num_instances){
2154 
2155  /*Allocate matrices :*/
2156  matrices=xNew<IssmDouble*>(num_instances);
2157  mdims=xNew<int>(num_instances);
2158  ndims=xNew<int>(num_instances);
2159 
2160  for(int i=0;i<num_instances;i++){
2161 
2162  if(my_rank==0){
2163  code=codes[i];
2164 
2165  if((code!=5) && (code!=6) && (code!=7))_error_("expecting a IssmDouble, integer or boolean matrix for \""<<data_name<<"\"");
2166 
2167  /*We have to read a matrix from disk. First read the dimensions of the matrix, then the whole matrix: */
2168  /*numberofelements: */
2169  fsetpos(fid,file_positions+i);
2170  if(fread(&M,sizeof(int),1,fid)!=1) _error_("could not read number of rows for matrix ");
2171  }
2173 
2174  if(my_rank==0){
2175  if(fread(&N,sizeof(int),1,fid)!=1) _error_("could not read number of columns for matrix ");
2176  }
2178 
2179  /*Now allocate matrix: */
2180  if(M*N){
2181  pmatrix=xNew<IssmPDouble>(M*N);
2182 
2183  /*Read matrix on node 0, then broadcast: */
2184  if(my_rank==0){
2185  if(fread(pmatrix,M*N*sizeof(IssmPDouble),1,fid)!=1) _error_("could not read matrix ");
2186  }
2188 
2189  //if(this->independents[data_enum]){ FIXME
2190  // /*this data has already been checked out! So cancel all that we've done here, and return
2191  // * the data[data_enum] directly: */
2192  // matrix=this->data[data_enum];
2193  //}
2194  //else{
2195  matrix=xNew<IssmDouble>(M*N);
2196  for (int i=0;i<M*N;++i) matrix[i]=pmatrix[i];
2197  //}
2198  xDelete<IssmPDouble>(pmatrix);
2199  }
2200  else
2201  matrix=NULL;
2202 
2203  /*Assign: */
2204  mdims[i]=M;
2205  matrices[i]=matrix;
2206  ndims[i]=N;
2207  }
2208  }
2209 
2210  /*Free ressources:*/
2211  xDelete<fpos_t>(file_positions);
2212  xDelete<int>(codes);
2213 
2214  /*Assign output pointers: */
2215  *pmatrices=matrices;
2216  if(pmdims){
2217  *pmdims=mdims;
2218  }
2219  else{
2220  xDelete<int>(mdims);
2221  }
2222  if(pndims){
2223  *pndims=ndims;
2224  }
2225  else{
2226  xDelete<int>(ndims);
2227  }
2228  *pnumrecords=num_instances;
2229 }

◆ FetchMultipleData() [3/5]

void IoModel::FetchMultipleData ( int ***  pmatrices,
int **  pmdims,
int **  pndims,
int *  pnumrecords,
const char *  data_name 
)

Definition at line 2231 of file IoModel.cpp.

2231  {/*{{{*/
2232 
2233  int num_instances;
2234  fpos_t* file_positions=NULL;
2235 
2236  /*output: */
2237  int **matrices = NULL;
2238  int *mdims = NULL;
2239  int *ndims = NULL;
2240 
2241  /*intermediary: */
2242  int M, N;
2243  IssmPDouble *pmatrix = NULL;
2244  IssmDouble *matrix = NULL;
2245  int *integer_matrix=NULL;
2246  int *codes = NULL;
2247  int code;
2248 
2249  /*recover my_rank:*/
2250  int my_rank=IssmComm::GetRank();
2251 
2252  /*Get file pointers to beginning of the data (multiple instances of it): */
2253  file_positions=this->SetFilePointersToData(&codes,NULL,&num_instances,data_name);
2254 
2255  if(num_instances){
2256 
2257  /*Allocate matrices :*/
2258  matrices=xNew<int*>(num_instances);
2259  mdims=xNew<int>(num_instances);
2260  ndims=xNew<int>(num_instances);
2261 
2262  for(int i=0;i<num_instances;i++){
2263 
2264  if(my_rank==0){
2265  code=codes[i];
2266 
2267  if((code!=5) && (code!=6) && (code!=7))_error_("expecting a IssmDouble, integer or boolean matrix for \""<<data_name<<"\"");
2268 
2269  /*We have to read a matrix from disk. First read the dimensions of the matrix, then the whole matrix: */
2270  /*numberofelements: */
2271  fsetpos(fid,file_positions+i);
2272  if(fread(&M,sizeof(int),1,fid)!=1) _error_("could not read number of rows for matrix ");
2273  }
2275 
2276  if(my_rank==0){
2277  if(fread(&N,sizeof(int),1,fid)!=1) _error_("could not read number of columns for matrix ");
2278  }
2280 
2281  /*Now allocate matrix: */
2282  if(M*N){
2283  pmatrix=xNew<IssmPDouble>(M*N);
2284  integer_matrix=xNew<int>(M*N);
2285 
2286  /*Read matrix on node 0, then broadcast: */
2287  if(my_rank==0){
2288  if(fread(pmatrix,M*N*sizeof(IssmPDouble),1,fid)!=1) _error_("could not read matrix ");
2289  }
2291 
2292  //if(this->independents[data_enum]){ FIXME
2293  // /*this data has already been checked out! So cancel all that we've done here, and return
2294  // * the data[data_enum] directly: */
2295  // matrix=this->data[data_enum];
2296  // for (int i=0;i<M*N;++i) integer_matrix[i]=reCast<int>(matrix[i]);
2297  //}
2298  //else{
2299  for (int i=0;i<M*N;++i) integer_matrix[i]=pmatrix[i];
2300  //}
2301  xDelete<IssmPDouble>(pmatrix);
2302  }
2303  else
2304  integer_matrix=NULL;
2305 
2306  /*Assign: */
2307  mdims[i]=M;
2308  matrices[i]=integer_matrix;
2309  ndims[i]=N;
2310  }
2311  }
2312 
2313  /*Free ressources:*/
2314  xDelete<fpos_t>(file_positions);
2315  xDelete<int>(codes);
2316 
2317  /*Assign output pointers: */
2318  *pmatrices=matrices;
2319  if(pmdims){
2320  *pmdims=mdims;
2321  }
2322  else{
2323  xDelete<int>(mdims);
2324  }
2325  if(pndims){
2326  *pndims=ndims;
2327  }
2328  else{
2329  xDelete<int>(ndims);
2330  }
2331  *pnumrecords=num_instances;
2332 }

◆ FetchMultipleData() [4/5]

void IoModel::FetchMultipleData ( int **  pvector,
int *  pnum_instances,
const char *  data_name 
)

Definition at line 2024 of file IoModel.cpp.

2024  {/*{{{*/
2025 
2026  int num_instances;
2027  fpos_t* file_positions=NULL;
2028 
2029  /*output: */
2030  int* vector=NULL;
2031 
2032  /*intermediary: */
2033  int integer;
2034  int *codes = NULL;
2035  int code;
2036 
2037  /*recover my_rank:*/
2038  int my_rank=IssmComm::GetRank();
2039 
2040  /*Get file pointers to beginning of the data (multiple instances of it): */
2041  file_positions=this->SetFilePointersToData(&codes,NULL,&num_instances,data_name);
2042 
2043  if(num_instances){
2044 
2045  /*Allocate vector :*/
2046  vector=xNew<int>(num_instances);
2047 
2048  for(int i=0;i<num_instances;i++){
2049 
2050  if(my_rank==0){
2051  code=codes[i];
2052 
2053  if(code!=2)_error_("expecting an integer for \""<<data_name<<"\"");
2054 
2055  /*We have to read a integer from disk. First read the dimensions of the integer, then the integer: */
2056  fsetpos(fid,file_positions+i);
2057  if(my_rank==0){
2058  if(fread(&integer,sizeof(int),1,fid)!=1) _error_("could not read integer ");
2059  }
2060  }
2062 
2063  /*Assign: */
2064  vector[i]=integer;
2065  }
2066  }
2067 
2068  /*Free ressources:*/
2069  xDelete<fpos_t>(file_positions);
2070  xDelete<int>(codes);
2071 
2072  /*Assign output pointers: */
2073  *pvector=vector;
2074  *pnum_instances=num_instances;
2075 }

◆ FetchMultipleData() [5/5]

void IoModel::FetchMultipleData ( IssmDouble **  pvector,
int *  pnum_instances,
const char *  data_name 
)

Definition at line 2077 of file IoModel.cpp.

2077  {/*{{{*/
2078 
2079  int num_instances;
2080  fpos_t* file_positions=NULL;
2081 
2082  /*output: */
2083  IssmDouble* vector=NULL;
2084 
2085  /*intermediary: */
2086  IssmPDouble scalar;
2087  int *codes = NULL;
2088  int code;
2089 
2090  /*recover my_rank:*/
2091  int my_rank=IssmComm::GetRank();
2092 
2093  /*Get file pointers to beginning of the data (multiple instances of it): */
2094  file_positions=this->SetFilePointersToData(&codes,NULL,&num_instances,data_name);
2095 
2096  if(num_instances){
2097 
2098  /*Allocate vector :*/
2099  vector=xNew<IssmDouble>(num_instances);
2100 
2101  for(int i=0;i<num_instances;i++){
2102 
2103  if(my_rank==0){
2104  code=codes[i];
2105 
2106  if(code!=3)_error_("expecting a double for \""<<data_name<<"\"");
2107 
2108  /*We have to read a double from disk: */
2109  fsetpos(fid,file_positions+i);
2110  if(my_rank==0){
2111  if(fread(&scalar,sizeof(IssmPDouble),1,fid)!=1) _error_("could not read scalar ");
2112  }
2113  }
2115 
2116  /*Assign: */
2117  vector[i]=scalar;
2118  }
2119  }
2120 
2121  /*Free ressources:*/
2122  xDelete<fpos_t>(file_positions);
2123  xDelete<int>(codes);
2124 
2125  /*Assign output pointers: */
2126  *pvector=vector;
2127  *pnum_instances=num_instances;
2128 }

◆ SetFilePointersToData()

fpos_t * IoModel::SetFilePointersToData ( int **  pcodes,
int **  pvector_types,
int *  pnum_instances,
const char *  data_name 
)

Definition at line 2471 of file IoModel.cpp.

2471  {/*{{{*/
2472 
2473  int found = 0;
2474  const char* mddot = "md.";
2475  char* record_name = NULL;
2476  int record_name_size;
2477  long long record_length;
2478  int record_code; //1 to 7 number
2479  int vector_type; //1 to 7 number
2480  int *vector_types = NULL;
2481  int *codes = NULL;
2482  int num_instances = 0;
2483  int counter;
2484  fpos_t *file_positions = NULL;
2485 
2486  /*recover my_rank:*/
2487  int my_rank=IssmComm::GetRank();
2488  _assert_(strncmp(data_name,mddot,3)==0);
2489 
2490  /*Go find in the binary file, the data we want to fetch and count the number of
2491  * instances it appears: */
2492  if(my_rank==0){
2493 
2494  /*First set FILE* position to the beginning of the file: */
2495  fseek(fid,0,SEEK_SET);
2496 
2497  /*Now march through file looking for the correct data identifier: */
2498  for(;;){
2499  /*Read size of first string name: */
2500  if(fread(&record_name_size,sizeof(int),1,fid)==0){
2501  /*we have reached the end of the file. break: */
2502  xDelete<char>(record_name);
2503  break;
2504  }
2505  if(record_name_size<3 || record_name_size>80){
2506  _error_("error while looking in binary file. Found a string of size "<<record_name_size);
2507  }
2508 
2509  /*Allocate string of correct size: */
2510  record_name=xNew<char>(record_name_size+1);
2511  record_name[record_name_size]='\0';
2512 
2513  /*Read record_name: */
2514  if(fread(record_name,record_name_size*sizeof(char),1,fid)==0){
2515  break;
2516  }
2517  if(strncmp(record_name,mddot,3)!=0){
2518  _error_("error while reading binary file: record does not start with \"md.\": "<<record_name);
2519  }
2520 
2521  /*Is this the record sought for? : */
2522  if(strcmp(record_name,data_name)==0) num_instances++;
2523 
2524  /*Read the record length, and use it to skip the record: */
2525  if(fread(&record_length,sizeof(long long),1,fid)!=1) _error_("Could not read record_length");
2526  fseek(fid,record_length,SEEK_CUR);
2527  xDelete<char>(record_name);
2528  }
2529 
2530  /*Ok, initialize the number of file handles we are going to return: */
2531  if(num_instances){
2532  file_positions = xNew<fpos_t>(num_instances);
2533  codes = xNew<int>(num_instances);
2534  vector_types = xNew<int>(num_instances);
2535  }
2536 
2537  /*Reset FILE* position to the beginning of the file, and start again, this time saving the data information
2538  * as we find it: */
2539  counter=0;
2540  fseek(fid,0,SEEK_SET);
2541 
2542  for(;;){
2543  /*Read size of first string name: */
2544  if(fread(&record_name_size,sizeof(int),1,fid)==0){
2545  /*we have reached the end of the file. break: */
2546  break;
2547  }
2548  if(record_name_size<3 || record_name_size>80){
2549  _error_("error while looking in binary file. Found a string of size "<<record_name_size);
2550  }
2551 
2552  /*Allocate string of correct size: */
2553  record_name=xNew<char>(record_name_size+1);
2554  record_name[record_name_size]='\0';
2555 
2556  /*Read record_name: */
2557  if(fread(record_name,record_name_size*sizeof(char),1,fid)==0){
2558  /*we have reached the end of the file. break: */
2559  xDelete<char>(record_name);
2560  break;
2561  }
2562  if(strncmp(record_name,mddot,3)!=0){
2563  _error_("error while reading binary file: record does not start with \"md.\": "<<record_name);
2564  }
2565 
2566  /*Is this the record sought for? : */
2567  if(strcmp(record_name,data_name)==0){
2568  /*Ok, we have found the correct string. Pass the record length, and read data type code: */
2569  fseek(fid,sizeof(long long),SEEK_CUR);
2570  if(fread(&record_code,sizeof(int),1,fid)!=1) _error_("Could not read record_code");
2571 
2572  /*if record_code points to a vector, get its type (nodal or elementary): */
2573  if(5<=record_code && record_code<=7){
2574  if(fread(&vector_type,sizeof(int),1,fid)!=1) _error_("Could not read vector_type");
2575  }
2576  codes[counter] = record_code;
2577  vector_types[counter] = vector_type;
2578  fgetpos(fid,file_positions+counter);
2579 
2580  /*backup and skip over the record, as we have more work to do: */
2581  if(5<=record_code && record_code<=7) fseek(fid,-sizeof(int),SEEK_CUR); /*rewind for nodal or elementary type*/
2582  fseek(fid,-sizeof(int),SEEK_CUR);/*rewind for data code*/
2583  fseek(fid,-sizeof(long long),SEEK_CUR);/*rewind for record length*/
2584 
2585  /*increment counter: */
2586  counter++;
2587  }
2588 
2589  /*Read the record length, and use it to skip this record, as it has already been processed: */
2590  if(fread(&record_length,sizeof(long long),1,fid)!=1) _error_("Could not read record_length");
2591  /*skip: */
2592  fseek(fid,record_length,SEEK_CUR);
2593  xDelete<char>(record_name);
2594  }
2595  }
2596 
2597  /*Broadcast data: */
2598  ISSM_MPI_Bcast(&num_instances,1,ISSM_MPI_INT,0,IssmComm::GetComm());
2599 
2600  /*Assign output pointers:*/
2601  *pcodes = codes;
2602  *pnum_instances = num_instances;
2603  if(pvector_types){
2604  *pvector_types=vector_types;
2605  }
2606  else{
2607  xDelete<int>(vector_types);
2608  }
2609  return file_positions;
2610 }

◆ SetFilePointerToData()

FILE * IoModel::SetFilePointerToData ( int *  pcode,
int *  pvector_type,
const char *  data_name 
)

Definition at line 2612 of file IoModel.cpp.

2612  {/*{{{*/
2613 
2614  int my_rank;
2615 
2616  int found = 0;
2617  const char* mddot = "md.";
2618  char* record_name = NULL;
2619  int record_name_size;
2620  long long record_length;
2621  int record_code; //1 to 7 number
2622  int vector_type = 0; //nodal or elementary
2623 
2624  /*recover my_rank:*/
2625  my_rank=IssmComm::GetRank();
2626  if(strncmp(data_name,mddot,3)!=0){
2627  _error_("Cannot fetch \""<<data_name<<"\" does not start with \""<<mddot<<"\"");
2628  }
2629 
2630  /*Go find in the binary file, the position of the data we want to fetch: */
2631  if(my_rank==0){
2632 
2633  /*First set FILE* position to the beginning of the file: */
2634  _assert_(fid);
2635  fseek(fid,0,SEEK_SET);
2636 
2637  /*Now march through file looking for the correct data identifier: */
2638  for(;;){
2639  /*Read size of first string name: */
2640  if(fread(&record_name_size,sizeof(int),1,fid)==0){
2641  /*we have reached the end of the file. break: */
2642  xDelete<char>(record_name);
2643  break;
2644  }
2645  if(record_name_size<3 || record_name_size>80){
2646  _error_("error while looking in binary file. Found a string of size "<<record_name_size);
2647  }
2648 
2649  /*Allocate string of correct size: */
2650  record_name=xNew<char>(record_name_size+1);
2651  record_name[record_name_size]='\0';
2652 
2653  /*Read record_name: */
2654  if(fread(record_name,record_name_size*sizeof(char),1,fid)==0){
2655  /*we have reached the end of the file. break: */
2656  found=0;
2657  xDelete<char>(record_name);
2658  break;
2659  }
2660  if(strncmp(record_name,mddot,3)!=0){
2661  _error_("error while reading binary file: record does not start with \"md.\": "<<record_name);
2662  }
2663 
2664  /*Is this the record sought for? : */
2665  if(strcmp(record_name,data_name)==0){
2666  /*Ok, we have found the correct string. Pass the record length, and read data type code: */
2667  fseek(fid,sizeof(long long),SEEK_CUR);
2668  if(fread(&record_code,sizeof(int),1,fid)!=1) _error_("Could not read record_code");
2669  /*if record_code points to a vector, get its type (nodal or elementary): */
2670  if((5<=record_code && record_code<=7) || record_code==10){
2671  if(fread(&vector_type,sizeof(int),1,fid)!=1) _error_("Could not read vector_type");
2672  }
2673  found=1;
2674  xDelete<char>(record_name);
2675  break;
2676  }
2677  else{
2678  /*This is not the correct string, read the record length, and use it to skip this record: */
2679  if(fread(&record_length,sizeof(long long),1,fid)!=1) _error_("Could not read record_length");
2680  /*skip: */
2681  fseek(fid,record_length,SEEK_CUR);
2682  xDelete<char>(record_name);
2683  }
2684  }
2685  }
2687  if(!found) _error_("could not find data with name \"" << data_name << "\" in binary file");
2688 
2689  /*Broadcast code and vector type: */
2690  ISSM_MPI_Bcast(&record_code,1,ISSM_MPI_INT,0,IssmComm::GetComm());
2691  ISSM_MPI_Bcast(&vector_type,1,ISSM_MPI_INT,0,IssmComm::GetComm());
2692 
2693  /*Assign output pointers:*/
2694  *pcode=record_code;
2695  if(pvector_type)*pvector_type=vector_type;
2696 
2697  return fid;
2698 }

◆ StartTrace()

void IoModel::StartTrace ( bool  trace)

Definition at line 2700 of file IoModel.cpp.

2700  {/*{{{*/
2701 
2702  bool autodiff = false;
2703  bool iscontrol = false;
2704  bool keep=false;
2705  IssmDouble gcTriggerRatio;
2706  IssmDouble gcTriggerMaxSize;
2707  IssmDouble obufsize;
2708  IssmDouble lbufsize;
2709  IssmDouble cbufsize;
2710  IssmDouble tbufsize;
2711 
2712  int my_rank=IssmComm::GetRank();
2713 
2714  this->FetchData(&autodiff,"md.autodiff.isautodiff");
2715  this->FetchData(&iscontrol,"md.inversion.iscontrol");
2716 
2717  if(trace || (autodiff && !iscontrol)){
2718 
2719  #if defined(_HAVE_ADOLC_)
2720  /*Retrieve parameters: */
2721  this->FetchData(&keep,"md.autodiff.keep");
2722  int keepTaylors=keep?1:0;
2723  this->FetchData(&gcTriggerRatio,"md.autodiff.gcTriggerRatio");
2724  this->FetchData(&gcTriggerMaxSize,"md.autodiff.gcTriggerMaxSize");
2725  this->FetchData(&obufsize,"md.autodiff.obufsize");
2726  this->FetchData(&lbufsize,"md.autodiff.lbufsize");
2727  this->FetchData(&cbufsize,"md.autodiff.cbufsize");
2728  this->FetchData(&tbufsize,"md.autodiff.tbufsize");
2729 
2730  /*Set garbage collection parameters: */
2731  setStoreManagerControl(reCast<IssmPDouble>(gcTriggerRatio),reCast<size_t>(gcTriggerMaxSize));
2732 
2733  /*Start trace: */
2734  int skipFileDeletion=1;
2735  trace_on(my_rank,keepTaylors,reCast<size_t>(obufsize),reCast<size_t>(lbufsize),reCast<size_t>(cbufsize),reCast<size_t>(tbufsize),skipFileDeletion);
2736 
2737  #elif defined(_HAVE_CODIPACK_)
2738  //fprintf(stderr, "*** Codipack IoModel::StartTrace\n");
2739  /*
2740  * FIXME codi
2741  * - ADOL-C variant uses fine grained tracing with various arguments
2742  * - ADOL-C variant sets a garbage collection parameter for its tape
2743  * -> These parameters are not read for the CoDiPack ISSM version!
2744  */
2745  auto& tape_codi = IssmDouble::getGlobalTape();
2746  tape_codi.setActive();
2747  #if _AD_TAPE_ALLOC_
2748  //alloc_profiler.Tag(StartInit, true);
2749  IssmDouble x_t(1.0), y_t(1.0);
2750  tape_codi.registerInput(y_t);
2751  int codi_allocn = 0;
2752  this->FetchData(&codi_allocn,"md.autodiff.tapeAlloc");
2753  for(int i = 0;i < codi_allocn;++i) {
2754  x_t = y_t * y_t;
2755  }
2756  /*
2757  std::stringstream out_s;
2758  IssmDouble::getGlobalTape().printStatistics(out_s);
2759  _printf0_("StartTrace::Tape Statistics : TapeAlloc count=[" << codi_allocn << "]\n" << out_s.str());
2760  */
2761  tape_codi.reset();
2762  //alloc_profiler.Tag(FinishInit, true);
2763  #endif
2764  #endif
2765  }
2766 
2767 }

Field Documentation

◆ constants

std::vector<IoConstant*> IoModel::constants
private

Definition at line 51 of file IoModel.h.

◆ data

std::vector<IoData*> IoModel::data
private

Definition at line 52 of file IoModel.h.

◆ fid

FILE* IoModel::fid

Definition at line 60 of file IoModel.h.

◆ solution_enum

int IoModel::solution_enum

Definition at line 63 of file IoModel.h.

◆ my_elements

bool* IoModel::my_elements

Definition at line 66 of file IoModel.h.

◆ my_faces

bool* IoModel::my_faces

Definition at line 67 of file IoModel.h.

◆ my_vfaces

bool* IoModel::my_vfaces

Definition at line 68 of file IoModel.h.

◆ my_edges

bool* IoModel::my_edges

Definition at line 69 of file IoModel.h.

◆ my_vedges

bool* IoModel::my_vedges

Definition at line 70 of file IoModel.h.

◆ my_hedges

bool* IoModel::my_hedges

Definition at line 71 of file IoModel.h.

◆ my_vertices

bool* IoModel::my_vertices

Definition at line 72 of file IoModel.h.

◆ my_vertices_lids

int* IoModel::my_vertices_lids

Definition at line 73 of file IoModel.h.

◆ epart

int* IoModel::epart

Definition at line 74 of file IoModel.h.

◆ domaindim

int IoModel::domaindim

Definition at line 77 of file IoModel.h.

◆ domaintype

int IoModel::domaintype

Definition at line 78 of file IoModel.h.

◆ elements

int* IoModel::elements

Definition at line 79 of file IoModel.h.

◆ edges

int* IoModel::edges

Definition at line 80 of file IoModel.h.

◆ verticaledges

int* IoModel::verticaledges

Definition at line 81 of file IoModel.h.

◆ horizontaledges

int* IoModel::horizontaledges

Definition at line 82 of file IoModel.h.

◆ elementtoedgeconnectivity

int* IoModel::elementtoedgeconnectivity

Definition at line 83 of file IoModel.h.

◆ elementtoverticaledgeconnectivity

int* IoModel::elementtoverticaledgeconnectivity

Definition at line 84 of file IoModel.h.

◆ elementtohorizontaledgeconnectivity

int* IoModel::elementtohorizontaledgeconnectivity

Definition at line 85 of file IoModel.h.

◆ elementtofaceconnectivity

int* IoModel::elementtofaceconnectivity

Definition at line 86 of file IoModel.h.

◆ elementtoverticalfaceconnectivity

int* IoModel::elementtoverticalfaceconnectivity

Definition at line 87 of file IoModel.h.

◆ faces

int* IoModel::faces

Definition at line 88 of file IoModel.h.

◆ verticalfaces

int* IoModel::verticalfaces

Definition at line 89 of file IoModel.h.

◆ facescols

int IoModel::facescols

Definition at line 90 of file IoModel.h.

◆ meshelementtype

int IoModel::meshelementtype

Definition at line 91 of file IoModel.h.

◆ numbernodetoelementconnectivity

int* IoModel::numbernodetoelementconnectivity

Definition at line 92 of file IoModel.h.

◆ numberofedges

int IoModel::numberofedges

Definition at line 93 of file IoModel.h.

◆ numberofverticaledges

int IoModel::numberofverticaledges

Definition at line 94 of file IoModel.h.

◆ numberofhorizontaledges

int IoModel::numberofhorizontaledges

Definition at line 95 of file IoModel.h.

◆ numberofelements

int IoModel::numberofelements

Definition at line 96 of file IoModel.h.

◆ numberoffaces

int IoModel::numberoffaces

Definition at line 97 of file IoModel.h.

◆ numberofverticalfaces

int IoModel::numberofverticalfaces

Definition at line 98 of file IoModel.h.

◆ numberofvertices

int IoModel::numberofvertices

Definition at line 99 of file IoModel.h.

◆ singlenodetoelementconnectivity

int* IoModel::singlenodetoelementconnectivity

Definition at line 100 of file IoModel.h.


The documentation for this class was generated from the following files:
GenericOption::size
int size[2]
Definition: GenericOption.h:27
IoModel::solution_enum
int solution_enum
Definition: IoModel.h:63
IoModel::AddConstantIndependent
void AddConstantIndependent(IoConstant *constant_in)
Definition: IoModel.cpp:305
Element::SetIntInput
void SetIntInput(Inputs2 *inputs2, int enum_in, int value)
Definition: Element.cpp:3362
IoData::isindependent
bool isindependent
Definition: IoModel.h:38
_assert_
#define _assert_(ignore)
Definition: exceptions.h:37
IoModel::verticalfaces
int * verticalfaces
Definition: IoModel.h:89
IssmDouble
double IssmDouble
Definition: types.h:37
Param
Definition: Param.h:21
_printf0_
#define _printf0_(StreamArgs)
Definition: Print.h:29
Element::SetBoolInput
void SetBoolInput(Inputs2 *inputs2, int enum_in, bool value)
Definition: Element.cpp:3355
IoModel::numbernodetoelementconnectivity
int * numbernodetoelementconnectivity
Definition: IoModel.h:92
IoModel::elementtohorizontaledgeconnectivity
int * elementtohorizontaledgeconnectivity
Definition: IoModel.h:85
IoModel::elementtofaceconnectivity
int * elementtofaceconnectivity
Definition: IoModel.h:86
IoModel::DeclareIndependents
void DeclareIndependents(bool trace, IssmPDouble *X)
Definition: IoModel.cpp:450
GenericOption::name
char * name
Definition: GenericOption.h:25
IoModel::my_vfaces
bool * my_vfaces
Definition: IoModel.h:68
IoModel::constants
std::vector< IoConstant * > constants
Definition: IoModel.h:51
IoModel::my_hedges
bool * my_hedges
Definition: IoModel.h:71
IssmComm::GetComm
static ISSM_MPI_Comm GetComm(void)
Definition: IssmComm.cpp:30
IoModel::my_vertices_lids
int * my_vertices_lids
Definition: IoModel.h:73
IoCodeToEnumMaterials
int IoCodeToEnumMaterials(int enum_in)
Definition: IoCodeConversions.cpp:258
IoModel::numberoffaces
int numberoffaces
Definition: IoModel.h:97
IoModel::my_elements
bool * my_elements
Definition: IoModel.h:66
Param::SetEnum
virtual void SetEnum(int enum_in)=0
IoModel::data
std::vector< IoData * > data
Definition: IoModel.h:52
Param::copy
virtual Param * copy()=0
GenericOption::value
OptionType value
Definition: GenericOption.h:26
IoConstant::constant
Param * constant
Definition: IoModel.h:21
IoModel::my_vertices
bool * my_vertices
Definition: IoModel.h:72
IoCodeToEnumHydrology
int IoCodeToEnumHydrology(int enum_in)
Definition: IoCodeConversions.cpp:248
IoData::data
IssmDouble * data
Definition: IoModel.h:37
Element
Definition: Element.h:41
IoModel::numberofvertices
int numberofvertices
Definition: IoModel.h:99
ISSM_MPI_DOUBLE
#define ISSM_MPI_DOUBLE
Definition: issmmpi.h:125
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::AddConstant
void AddConstant(IoConstant *constant_in)
Definition: IoModel.cpp:288
ISSM_MPI_INT
#define ISSM_MPI_INT
Definition: issmmpi.h:127
IoModel::numberofelements
int numberofelements
Definition: IoModel.h:96
IoModel::facescols
int facescols
Definition: IoModel.h:90
IoData
Definition: IoModel.h:34
IoCodeToEnumNature
int IoCodeToEnumNature(int enum_in)
Definition: IoCodeConversions.cpp:270
IoModel::elementtoverticaledgeconnectivity
int * elementtoverticaledgeconnectivity
Definition: IoModel.h:84
IoModel::AddDataIndependent
void AddDataIndependent(IoData *data_in)
Definition: IoModel.cpp:333
IoModel::verticaledges
int * verticaledges
Definition: IoModel.h:81
IoModel::elementtoedgeconnectivity
int * elementtoedgeconnectivity
Definition: IoModel.h:83
IoModel::FindConstant
void FindConstant(bool *pvalue, const char *constant_name)
Definition: IoModel.cpp:2362
IoModel::horizontaledges
int * horizontaledges
Definition: IoModel.h:82
IoModel::my_edges
bool * my_edges
Definition: IoModel.h:69
IoModel::fid
FILE * fid
Definition: IoModel.h:60
IoModel::SetFilePointersToData
fpos_t * SetFilePointersToData(int **pcodes, int **pvector_types, int *pnum_instances, const char *data_name)
Definition: IoModel.cpp:2471
IoCodeToEnumSMB
int IoCodeToEnumSMB(int enum_in)
Definition: IoCodeConversions.cpp:199
IoModel::FetchData
void FetchData(bool *pboolean, const char *data_name)
Definition: IoModel.cpp:933
IoModel::numberofhorizontaledges
int numberofhorizontaledges
Definition: IoModel.h:95
IoModel::domaindim
int domaindim
Definition: IoModel.h:77
ISSM_MPI_PDOUBLE
#define ISSM_MPI_PDOUBLE
Definition: issmmpi.h:126
IoModel::FetchIndependentConstant
void FetchIndependentConstant(int *pXcount, IssmPDouble *X, const char *name)
Definition: IoModel.cpp:1814
ISSM_MPI_Bcast
int ISSM_MPI_Bcast(void *buffer, int count, ISSM_MPI_Datatype datatype, int root, ISSM_MPI_Comm comm)
Definition: issmmpi.cpp:162
IoConstant
Definition: IoModel.h:19
IoData::N
int N
Definition: IoModel.h:40
StringToEnumx
int StringToEnumx(const char *string_in, bool notfounderror=true)
Definition: StringToEnumx.cpp:14
IoModel::singlenodetoelementconnectivity
int * singlenodetoelementconnectivity
Definition: IoModel.h:100
Element::SetElementInput
virtual void SetElementInput(int enum_in, IssmDouble values)
Definition: Element.h:333
IoModel::CheckFile
void CheckFile(void)
Definition: IoModel.cpp:344
IoModel::FetchIndependentData
void FetchIndependentData(int *pXcount, IssmPDouble *X, const char *name)
Definition: IoModel.cpp:1868
IoData::name
char * name
Definition: IoModel.h:41
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
IssmComm::GetRank
static int GetRank(void)
Definition: IssmComm.cpp:34
ISSM_MPI_CHAR
#define ISSM_MPI_CHAR
Definition: issmmpi.h:124
_error_
#define _error_(StreamArgs)
Definition: exceptions.h:49
IoCodeToEnumFrontalforcings
int IoCodeToEnumFrontalforcings(int enum_in)
Definition: IoCodeConversions.cpp:241
IoModel::StartTrace
void StartTrace(bool trace)
Definition: IoModel.cpp:2700
IoModel::faces
int * faces
Definition: IoModel.h:88
IoModel::numberofverticalfaces
int numberofverticalfaces
Definition: IoModel.h:98
IoModel::numberofverticaledges
int numberofverticaledges
Definition: IoModel.h:94
IoModel::NumIndependents
int NumIndependents()
Definition: IoModel.cpp:2448
IoModel::SetFilePointerToData
FILE * SetFilePointerToData(int *pcode, int *pvector_type, const char *data_name)
Definition: IoModel.cpp:2612
IoData::code
int code
Definition: IoModel.h:36
IoData::layout
int layout
Definition: IoModel.h:39
IoModel::AddData
void AddData(IoData *data_in)
Definition: IoModel.cpp:316
IoModel::epart
int * epart
Definition: IoModel.h:74
IoModel::elements
int * elements
Definition: IoModel.h:79
IoModel::edges
int * edges
Definition: IoModel.h:80
IoModel::meshelementtype
int meshelementtype
Definition: IoModel.h:91
xMemCpy
T * xMemCpy(T *dest, const T *src, unsigned int size)
Definition: MemOps.h:152
IoModel::elementtoverticalfaceconnectivity
int * elementtoverticalfaceconnectivity
Definition: IoModel.h:87
IssmPDouble
IssmDouble IssmPDouble
Definition: types.h:38
IoModel::my_faces
bool * my_faces
Definition: IoModel.h:67
IoModel::domaintype
int domaintype
Definition: IoModel.h:78
IoConstant::name
char * name
Definition: IoModel.h:23
IoCodeToEnumTimestepping
int IoCodeToEnumTimestepping(int enum_in)
Definition: IoCodeConversions.cpp:282
Options::AddOption
int AddOption(Option *in_oobject)
Definition: Options.cpp:33
ISSM_MPI_LONG_LONG_INT
#define ISSM_MPI_LONG_LONG_INT
Definition: issmmpi.h:128
IoData::M
int M
Definition: IoModel.h:40
IoModel::FetchConstants
void FetchConstants(void)
Definition: IoModel.cpp:557
IoCodeToEnumAmr
int IoCodeToEnumAmr(int enum_in)
Definition: IoCodeConversions.cpp:289
IoCodeToEnumBasal
int IoCodeToEnumBasal(int enum_in)
Definition: IoCodeConversions.cpp:216
IoModel::numberofedges
int numberofedges
Definition: IoModel.h:93
GenericOption
Definition: GenericOption.h:22
Param::GetParameterValue
virtual void GetParameterValue(bool *pbool)=0
IoModel::my_vedges
bool * my_vedges
Definition: IoModel.h:70
IoCodeToEnumCalving
int IoCodeToEnumCalving(int enum_in)
Definition: IoCodeConversions.cpp:229