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

Last change on this file since 3775 was 3775, checked in by Eric.Larour, 15 years ago

Created include.h header file

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