source: issm/trunk/src/c/parallel/transient.cpp@ 3086

Last change on this file since 3086 was 3086, checked in by seroussi, 15 years ago

added slopecompute solution in parallel

File size: 5.3 KB
Line 
1/*!\file: transient.cpp
2 * \brief: transient solution
3 */
4
5#include "../issm.h"
6#include "./parallel.h"
7
8#undef __FUNCT__
9#define __FUNCT__ "transient"
10
11#ifdef HAVE_CONFIG_H
12 #include "config.h"
13#else
14#error "Cannot compile with HAVE_CONFIG_H symbol! run configure first!"
15#endif
16
17
18int main(int argc,char* *argv){
19
20 /*I/O: */
21 FILE* fid=NULL;
22 char* inputfilename=NULL;
23 char* outputfilename=NULL;
24 char* lockname=NULL;
25 int numberofnodes;
26 int qmu_analysis=0;
27 double waitonlock=0;
28
29 /*Model: */
30 Model* model=NULL;
31 int dim=-1;
32
33 /*Results: */
34 DataSet* results=NULL;
35 DataSet* processed_results=NULL;
36 Result* result=NULL;
37
38 ParameterInputs* inputs=NULL;
39
40 /*inputs: */
41 double* u_g=NULL;
42 double* m_g=NULL;
43 double* a_g=NULL;
44 double dt;
45 double yts;
46 Param* param=NULL;
47
48 /*time*/
49 double start, finish;
50 double start_core, finish_core;
51 double start_init, finish_init;
52
53 MODULEBOOT();
54
55 #if !defined(_PARALLEL_) || (defined(_PARALLEL_) && !defined(_HAVE_PETSC_))
56 throw ErrorException(__FUNCT__," parallel executable was compiled without support of parallel libraries!");
57 #endif
58
59 /*Initialize Petsc and get start time*/
60 PetscInitialize(&argc,&argv,(char *)0,"");
61 MPI_Barrier(MPI_COMM_WORLD); start=MPI_Wtime();
62
63 /*Size and rank: */
64 MPI_Comm_rank(MPI_COMM_WORLD,&my_rank);
65 MPI_Comm_size(MPI_COMM_WORLD,&num_procs);
66
67 _printf_("recover , input file name and output file name:\n");
68 inputfilename=argv[2];
69 outputfilename=argv[3];
70 lockname=argv[4];
71
72 /*Open handle to data on disk: */
73 fid=pfopen(inputfilename,"rb");
74
75 /*Initialize model structure: */
76 MPI_Barrier(MPI_COMM_WORLD); start_init=MPI_Wtime();
77 model=new Model();
78
79 _printf_("read and create finite element model:\n");
80 _printf_("\n reading diagnostic horiz model data:\n");
81 model->AddFormulation(fid,DiagnosticAnalysisEnum(),HorizAnalysisEnum());
82
83 _printf_("\n reading diagnostic vert model data:\n");
84 model->AddFormulation(fid,DiagnosticAnalysisEnum(),VertAnalysisEnum());
85
86 _printf_("\n reading diagnostic stokes model data:\n");
87 model->AddFormulation(fid,DiagnosticAnalysisEnum(),StokesAnalysisEnum());
88
89 _printf_("\n reading diagnostic hutter model data:\n");
90 model->AddFormulation(fid,DiagnosticAnalysisEnum(),HutterAnalysisEnum());
91
92 _printf_("\n reading surface and bed slope computation model data:\n");
93 model->AddFormulation(fid,SlopecomputeAnalysisEnum());
94
95 _printf_("\n reading prognositc model data:\n");
96 model->AddFormulation(fid,PrognosticAnalysisEnum());
97
98 /*Do we run in 3d?, in which case we need thermal and melting also:*/
99 model->FindParam(&dim,"dim");
100 if(dim==3){
101 _printf_("read and create thermal finite element model:\n");
102 model->AddFormulation(fid,ThermalAnalysisEnum(),TransientAnalysisEnum());
103 _printf_("read and create melting finite element model:\n");
104 model->AddFormulation(fid,MeltingAnalysisEnum(),TransientAnalysisEnum());
105 }
106
107 /*recover parameters: */
108 model->FindParam(&waitonlock,"waitonlock");
109 model->FindParam(&qmu_analysis,"qmu_analysis");
110
111 _printf_("initialize inputs:\n");
112 model->FindParam(&u_g,NULL,NULL,"u_g",PrognosticAnalysisEnum());
113 model->FindParam(&m_g,NULL,NULL,"m_g",PrognosticAnalysisEnum());
114 model->FindParam(&a_g,NULL,NULL,"a_g",PrognosticAnalysisEnum());
115 model->FindParam(&numberofnodes,"numberofnodes");
116 model->FindParam(&dt,"dt");
117 model->FindParam(&yts,"yts");
118
119
120 inputs=new ParameterInputs;
121 inputs->Add("velocity",u_g,3,numberofnodes);
122 inputs->Add("melting",m_g,1,numberofnodes);
123 inputs->Add("accumulation",a_g,1,numberofnodes);
124 inputs->Add("dt",dt*yts);
125
126 _printf_("initialize results:\n");
127 results=new DataSet(ResultsEnum());
128 MPI_Barrier(MPI_COMM_WORLD); finish_init=MPI_Wtime();
129
130 /*are we running the solution sequence, or a qmu wrapper around it? : */
131 if(!qmu_analysis){
132
133 /*run diagnostic analysis: */
134 _printf_("call computational core:\n");
135 MPI_Barrier(MPI_COMM_WORLD); start_core=MPI_Wtime( );
136 transient_core(results,model,inputs);
137 MPI_Barrier(MPI_COMM_WORLD); finish_core=MPI_Wtime( );
138
139 /*Add analysis_type to results: */
140 result=new Result(results->Size()+1,0,1,"analysis_type","transient");
141 results->AddObject(result);
142
143 _printf_("process results:\n");
144 ProcessResults(&processed_results,results,model,TransientAnalysisEnum());
145
146 _printf_("write results to disk:\n");
147 OutputResults(processed_results,outputfilename);
148 }
149 else{
150 /*run qmu analysis: */
151 _printf_("calling qmu analysis on transient core:\n");
152
153 #ifdef _HAVE_DAKOTA_
154 MPI_Barrier(MPI_COMM_WORLD); start_core=MPI_Wtime( );
155 Qmux(model,inputs,TransientAnalysisEnum(),NoneAnalysisEnum());
156 MPI_Barrier(MPI_COMM_WORLD); finish_core=MPI_Wtime( );
157 #else
158 throw ErrorException(__FUNCT__," Dakota not present, cannot do qmu!");
159 #endif
160 }
161
162 if (waitonlock>0){
163 _printf_("write lock file:\n");
164 WriteLockFile(lockname);
165 }
166
167 /*Free ressources:*/
168 delete results;
169 delete processed_results;
170 delete model;
171 delete inputs;
172 xfree((void**)&u_g);
173 xfree((void**)&m_g);
174 xfree((void**)&a_g);
175
176 /*Get finish time and close*/
177 MPI_Barrier(MPI_COMM_WORLD); finish = MPI_Wtime( );
178 _printf_("\n %-34s %f seconds \n","Model initialization elapsed time:",finish_init-start_init);
179 _printf_(" %-34s %f seconds \n","Core solution elapsed time:",finish_core-start_core);
180 _printf_(" %-34s %f seconds\n\n","Total elapsed time:",finish-start);
181 _printf_("closing MPI and Petsc\n");
182 PetscFinalize();
183
184 /*end module: */
185 MODULEEND();
186
187 return 0; //unix success return;
188}
Note: See TracBrowser for help on using the repository browser.