Ice Sheet System Model  4.18
Code documentation
issm_slr.cpp
Go to the documentation of this file.
1 
5 #include "./issm.h"
6 #include <stdlib.h>
7 
8 int main(int argc,char **argv){
9 
10  /*diverse:*/
11  int nummodels;
12  int* commsizes=NULL;
13  int* rankzeros=NULL;
14  char** dirnames=NULL;
15  char** modelnames=NULL;
16  int modelid;
17  int earthid;
18  int my_rank;
19  int count=0;
20  ISSM_MPI_Comm worldcomm;
21  ISSM_MPI_Comm modelcomm;
22  ISSM_MPI_Comm toearthcomm;
23  ISSM_MPI_Comm* fromicecomms=NULL;
24 
25  /*Initialize exception trapping: */
27 
28  /*Initialize environment (MPI, PETSC, MUMPS, etc ...)*/
29  worldcomm=EnvironmentInit(argc,argv);
30 
31  /*What is my rank?:*/
32  ISSM_MPI_Comm_rank(worldcomm,&my_rank);
33 
34  /*How many models are we going to run (along with description and number of dedicated cores):{{{*/
35  nummodels=(int) strtol(argv[4], (char **)NULL, 10);
36  commsizes=xNew<int>(nummodels);
37  dirnames=xNew<char*>(nummodels);
38  modelnames=xNew<char*>(nummodels);
39  rankzeros=xNew<int>(nummodels);
40  for(int i=0;i<nummodels;i++){
41  char* string=NULL;
42 
43  string=xNew<char>(strlen(argv[5+3*i])+1);
44  xMemCpy<char>(string,argv[5+3*i],strlen(argv[5+3*i])+1);
45  dirnames[i]=string;
46 
47  string=xNew<char>(strlen(argv[5+3*i+1])+1);
48  xMemCpy<char>(string,argv[5+3*i+1],strlen(argv[5+3*i+1])+1);
49  modelnames[i]=string;
50 
51  commsizes[i]=(int) strtol(argv[5+3*i+2], (char **)NULL, 10);
52  }
53 
54  /*Figure out which model each cpu will belong to: */
55  count=0;
56  for(int i=0;i<nummodels;i++){
57  if(my_rank>=count && my_rank<(count+commsizes[i])){
58  modelid=i;
59  break;
60  }
61  count+=commsizes[i];
62  }
63  /*Buil array of who is rank 0 of their own group:*/
64  count=0;
65  for(int i=0;i<nummodels;i++){
66  rankzeros[i]=count;
67  count+=commsizes[i];
68  }
69  /*}}}*/
70 
71  /*Split world into sub-communicators for each and every model:*/
72  ISSM_MPI_Comm_split(worldcomm,modelid, my_rank, &modelcomm);
73 
74  /*Build inter communicators:*/
75  earthid=nummodels-1; //last model to be provided in the argument list if the earth model.
76  if(modelid==earthid){
77  fromicecomms=xNew<ISSM_MPI_Comm>(nummodels-1);
78  for(int i=0;i<earthid;i++){
79  ISSM_MPI_Intercomm_create( modelcomm, 0, worldcomm, rankzeros[i], i, fromicecomms+i); //communicate from local erth comm 9rank 0) to ice comm (rank 0) using modelid tag.
80  }
81  }
82  else{
83  ISSM_MPI_Intercomm_create( modelcomm, 0, worldcomm, rankzeros[earthid], modelid, &toearthcomm); //communicate from local ice comm (rank 0) to earth comm (rank 0) using modelid tag.
84  }
85 
86  /*Supply specific argc and argv for each sub-communicator (corresponding to each model specificatiions):{{{*/
87  char** arguments=xNew<char*>(4);
88  arguments[0]=xNew<char>(strlen(argv[0])+1); xMemCpy<char>(arguments[0],argv[0],strlen(argv[0])+1); //executable name
89  arguments[1]=xNew<char>(strlen(argv[1])+1); xMemCpy<char>(arguments[1],argv[1],strlen(argv[1])+1); //solution name
90  arguments[2]=xNew<char>(strlen(argv[5+3*modelid])+1); xMemCpy<char>(arguments[2],argv[5+3*modelid],strlen(argv[5+3*modelid])+1); //directory name
91  arguments[3]=xNew<char>(strlen(argv[5+3*modelid+1])+1); xMemCpy<char>(arguments[3],argv[5+3*modelid+1],strlen(argv[5+3*modelid+1])+1); //model name
92  /*}}}*/
93 
94  /*Initialize femmodel from arguments provided command line: */
95  FemModel *femmodel = new FemModel(4,arguments,modelcomm);
96 
97  /*Now that the models are initialized, keep communicator information in the parameters datasets of each model: */
102  if(modelid==earthid) femmodel->parameters->AddObject(new GenericParam<ISSM_MPI_Comm*>(fromicecomms,IcecapToEarthCommEnum));
104 
105  /*Solve: */
106  femmodel->Solve();
107 
108  /*Output results: */
110 
111  /*Wrap up: */
112  femmodel->CleanUp();
113 
114  /*Delete Model: */
115  delete femmodel;
116 
117  /*Finalize environment:*/
119 
120  /*Finalize exception trapping: */
122 
123  /*Free ressources:*/
124  xDelete<int>(commsizes);
125  for(int i=0;i<nummodels;i++){
126  char* string=NULL;
127  string=dirnames[i]; xDelete<char>(string);
128  string=modelnames[i]; xDelete<char>(string);
129  }
130  xDelete<char*>(dirnames);
131  xDelete<char*>(modelnames);
132 
133  /*Return unix success: */
134  return 0;
135 }
ISSM_MPI_Comm_split
int ISSM_MPI_Comm_split(ISSM_MPI_Comm comm, int color, int key, ISSM_MPI_Comm *newcomm)
Definition: issmmpi.cpp:528
ExceptionTrapEnd
#define ExceptionTrapEnd()
Definition: exceptions.h:64
FemModel::Solve
void Solve(void)
Definition: FemModel.cpp:883
OutputResultsx
void OutputResultsx(FemModel *femmodel)
Definition: OutputResultsx.cpp:17
ExceptionTrapBegin
#define ExceptionTrapBegin()
Definition: exceptions.h:61
FemModel::parameters
Parameters * parameters
Definition: FemModel.h:46
EnvironmentInit
ISSM_MPI_Comm EnvironmentInit(int argc, char **argv)
Definition: EnvironmentInit.cpp:12
WorldCommEnum
@ WorldCommEnum
Definition: EnumDefinitions.h:459
IcecapToEarthCommEnum
@ IcecapToEarthCommEnum
Definition: EnumDefinitions.h:200
Parameters::AddObject
void AddObject(Param *newparam)
Definition: Parameters.cpp:67
ISSM_MPI_Comm_rank
int ISSM_MPI_Comm_rank(ISSM_MPI_Comm comm, int *rank)
Definition: issmmpi.cpp:198
ModelIdEnum
@ ModelIdEnum
Definition: EnumDefinitions.h:274
EarthIdEnum
@ EarthIdEnum
Definition: EnumDefinitions.h:129
IntParam
Definition: IntParam.h:20
ISSM_MPI_Intercomm_create
int ISSM_MPI_Intercomm_create(ISSM_MPI_Comm comm, int local_leader, ISSM_MPI_Comm peer_comm, int remote_leader, int tag, ISSM_MPI_Comm *newintercomm)
Definition: issmmpi.cpp:542
FemModel
Definition: FemModel.h:31
EnvironmentFinalize
void EnvironmentFinalize(void)
Definition: EnvironmentFinalize.cpp:12
main
int main(int argc, char **argv)
Definition: issm_slr.cpp:8
ISSM_MPI_Comm
int ISSM_MPI_Comm
Definition: issmmpi.h:118
FemModel::CleanUp
void CleanUp(void)
Definition: FemModel.cpp:279
issm.h
prototype wrapper for issm.h
NumModelsEnum
@ NumModelsEnum
Definition: EnumDefinitions.h:276
GenericParam
Definition: GenericParam.h:26
femmodel
FemModel * femmodel
Definition: esmfbinders.cpp:16