source: issm/trunk/src/c/solutions/diagnostic.cpp@ 3895

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

Moved parallel to solutions

File size: 4.4 KB
Line 
1/*!\file: diagnostic.cpp
2 * \brief: diagnostic 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 "../modules.h"
17#include "./solutions.h"
18
19int main(int argc,char* *argv){
20
21 /*I/O: */
22 FILE* fid=NULL;
23 char* inputfilename=NULL;
24 char* outputfilename=NULL;
25 char* lockname=NULL;
26 bool qmu_analysis=false;
27 bool control_analysis=false;
28
29 /*Model: */
30 Model* model=NULL;
31
32 /*Results: */
33 DataSet* results=NULL;
34 DataSet* processed_results=NULL;
35 Result* result=NULL;
36
37 bool waitonlock=false;
38
39 /*time*/
40 double start, finish;
41 double start_core, finish_core;
42 double start_init, finish_init;
43
44 MODULEBOOT();
45
46 #if !defined(_PARALLEL_) || (defined(_PARALLEL_) && !defined(_HAVE_PETSC_))
47 ISSMERROR(" parallel executable was compiled without support of parallel libraries!");
48 #endif
49
50 /*Initialize Petsc and get start time*/
51 PetscInitialize(&argc,&argv,(char *)0,"");
52 MPI_Barrier(MPI_COMM_WORLD); start=MPI_Wtime();
53
54 /*Size and rank: */
55 MPI_Comm_rank(MPI_COMM_WORLD,&my_rank);
56 MPI_Comm_size(MPI_COMM_WORLD,&num_procs);
57
58 _printf_("recover , input file name and output file name:\n");
59 inputfilename=argv[2];
60 outputfilename=argv[3];
61 lockname=argv[4];
62
63 /*Initialize model structure: */
64 MPI_Barrier(MPI_COMM_WORLD); start_init=MPI_Wtime();
65 model=new Model();
66
67 /*Open handle to data on disk: */
68 fid=pfopen(inputfilename,"rb");
69
70 _printf_("read and create finite element model:\n");
71 _printf_("\n reading diagnostic horiz model data:\n");
72 model->AddFormulation(fid,DiagnosticAnalysisEnum,HorizAnalysisEnum);
73
74 _printf_("\n reading diagnostic vert model data:\n");
75 model->AddFormulation(fid,DiagnosticAnalysisEnum,VertAnalysisEnum);
76
77 _printf_("\n reading diagnostic stokes model data:\n");
78 model->AddFormulation(fid,DiagnosticAnalysisEnum,StokesAnalysisEnum);
79
80 _printf_("\n reading diagnostic hutter model data:\n");
81 model->AddFormulation(fid,DiagnosticAnalysisEnum,HutterAnalysisEnum);
82
83 _printf_("\n reading surface and bed slope computation model data:\n");
84 model->AddFormulation(fid,SlopecomputeAnalysisEnum);
85
86 /*get parameters: */
87 model->FindParam(&qmu_analysis,QmuAnalysisEnum);
88 model->FindParam(&control_analysis,ControlAnalysisEnum);
89 model->FindParam(&waitonlock,WaitOnLockEnum);
90
91 MPI_Barrier(MPI_COMM_WORLD); finish_init=MPI_Wtime();
92
93 /*are we running the solution sequence, or a qmu wrapper around it? : */
94 if(!qmu_analysis){
95 if(!control_analysis){
96
97 _printf_("call computational core:\n");
98 MPI_Barrier(MPI_COMM_WORLD); start_core=MPI_Wtime( );
99 results=diagnostic_core(model);
100 MPI_Barrier(MPI_COMM_WORLD); finish_core=MPI_Wtime( );
101
102 _printf_("process results:\n");
103 ProcessResults(&processed_results,results,model,DiagnosticAnalysisEnum);
104 }
105 else{
106 /*run control analysis: */
107 _printf_("call computational core:\n");
108 MPI_Barrier(MPI_COMM_WORLD); start_core=MPI_Wtime( );
109 results=control_core(model);
110 MPI_Barrier(MPI_COMM_WORLD); finish_core=MPI_Wtime( );
111
112 _printf_("process results:\n");
113 ProcessResults(&processed_results,results,model,ControlAnalysisEnum);
114 }
115
116 _printf_("write results to disk:\n");
117 OutputResults(processed_results,outputfilename);
118 }
119 else{
120 /*run qmu analysis: */
121 _printf_("calling qmu analysis on diagnostic core:\n");
122
123 #ifdef _HAVE_DAKOTA_
124 MPI_Barrier(MPI_COMM_WORLD); start_core=MPI_Wtime( );
125 Qmux(model,DiagnosticAnalysisEnum,NoneAnalysisEnum);
126 MPI_Barrier(MPI_COMM_WORLD); finish_core=MPI_Wtime( );
127 #else
128 ISSMERROR(" Dakota not present, cannot do qmu!");
129 #endif
130 }
131
132 if (waitonlock>0){
133 _printf_("write lock file:\n");
134 WriteLockFile(lockname);
135 }
136
137 /*Free ressources */
138 delete model;
139 delete results;
140 delete processed_results;
141
142 /*Get finish time and close*/
143 MPI_Barrier(MPI_COMM_WORLD); finish = MPI_Wtime( );
144 _printf_("\n %-34s %f seconds \n","Model initialization elapsed time:",finish_init-start_init);
145 _printf_(" %-34s %f seconds \n","Core solution elapsed time:",finish_core-start_core);
146 _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);
147 _printf_("closing MPI and Petsc\n");
148 PetscFinalize();
149
150 /*end module: */
151 MODULEEND();
152
153 return 0; //unix success return;
154}
Note: See TracBrowser for help on using the repository browser.