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

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

fixed leak

File size: 4.4 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
28 /*Model: */
29 Model* model=NULL;
30 int dim=-1;
31
32 /*Results: */
33 DataSet* results=NULL;
34 Result* result=NULL;
35
36 ParameterInputs* inputs=NULL;
37 int waitonlock=0;
38
39 /*inputs: */
40 double* u_g=NULL;
41 double* m_g=NULL;
42 double* a_g=NULL;
43 double dt;
44 double yts;
45 Param* param=NULL;
46
47 MODULEBOOT();
48
49 #if !defined(_PARALLEL_) || (defined(_PARALLEL_) && !defined(_HAVE_PETSC_))
50 throw ErrorException(__FUNCT__," parallel executable was compiled without support of parallel libraries!");
51 #endif
52
53 PetscInitialize(&argc,&argv,(char *)0,"");
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 model=new Model();
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 _printf_("\n reading prognositc model data:\n");
87 model->AddFormulation(fid,PrognosticAnalysisEnum());
88
89 /*Do we run in 3d?, in which case we need thermal and melting also:*/
90 model->FindParam(&dim,"dim");
91 if(dim==3){
92 _printf_("read and create thermal finite element model:\n");
93 model->AddFormulation(fid,ThermalAnalysisEnum(),TransientAnalysisEnum());
94 _printf_("read and create melting finite element model:\n");
95 model->AddFormulation(fid,MeltingAnalysisEnum(),TransientAnalysisEnum());
96 }
97
98 /*recover parameters: */
99 model->FindParam(&waitonlock,"waitonlock");
100 model->FindParam(&qmu_analysis,"qmu_analysis");
101
102 _printf_("initialize inputs:\n");
103 model->FindParam(&u_g,"u_g",PrognosticAnalysisEnum());
104 model->FindParam(&m_g,"m_g",PrognosticAnalysisEnum());
105 model->FindParam(&a_g,"a_g",PrognosticAnalysisEnum());
106 model->FindParam(&numberofnodes,"numberofnodes");
107 model->FindParam(&dt,"dt");
108 model->FindParam(&yts,"yts");
109
110
111 inputs=new ParameterInputs;
112 inputs->Add("velocity",u_g,3,numberofnodes);
113 inputs->Add("melting",m_g,1,numberofnodes);
114 inputs->Add("accumulation",a_g,1,numberofnodes);
115 inputs->Add("dt",dt*yts);
116
117 _printf_("initialize results:\n");
118 results=new DataSet(ResultsEnum());
119
120 /*are we running the solution sequence, or a qmu wrapper around it? : */
121 if(!qmu_analysis){
122
123 /*run diagnostic analysis: */
124 _printf_("call computational core:\n");
125 transient_core(results,model,inputs);
126 }
127 else{
128 /*run qmu analysis: */
129 _printf_("calling qmu analysis on transient core:\n");
130
131 #ifdef _HAVE_DAKOTA_
132 Qmux(model,inputs,TransientAnalysisEnum(),NoneAnalysisEnum());
133 #else
134 throw ErrorException(__FUNCT__," Dakota not present, cannot do qmu!");
135 #endif
136 }
137
138 /*Add analysis_type to results: */
139 result=new Result(results->Size()+1,0,1,"analysis_type","transient");
140 results->AddObject(result);
141
142 _printf_("process results:\n");
143 ProcessResults(&results,model,TransientAnalysisEnum());
144
145 _printf_("write results to disk:\n");
146 OutputResults(results,outputfilename);
147
148 _printf_("write lock file:\n");
149 if (waitonlock){
150 WriteLockFile(lockname);
151 }
152 /*Free ressources:*/
153 delete results;
154 delete model;
155 delete inputs;
156 xfree((void**)&u_g);
157 xfree((void**)&m_g);
158 xfree((void**)&a_g);
159
160 _printf_("closing MPI and Petsc\n");
161 PetscFinalize();
162
163 /*end module: */
164 MODULEEND();
165
166 return 0; //unix success return;
167}
Note: See TracBrowser for help on using the repository browser.