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

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

Added control type and analysis type in results

File size: 4.1 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 /*Fem models : */
29 FemModel femmodels[8];
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 Param* param=NULL;
45
46 MODULEBOOT();
47
48 #if !defined(_PARALLEL_) || (defined(_PARALLEL_) && !defined(_HAVE_PETSC_))
49 throw ErrorException(__FUNCT__," parallel executable was compiled without support of parallel libraries!");
50 #endif
51
52 PetscInitialize(&argc,&argv,(char *)0,"");
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 /*Open handle to data on disk: */
64 fid=pfopen(inputfilename,"rb");
65
66 _printf_("read and create finite element model:\n");
67 _printf_("\n reading diagnostic horiz model data:\n");
68 CreateFemModel(&femmodels[0],fid,"diagnostic","horiz");
69 _printf_("\n reading diagnostic vert model data:\n");
70 CreateFemModel(&femmodels[1],fid,"diagnostic","vert");
71 _printf_("\n reading diagnostic stokes model data:\n");
72 CreateFemModel(&femmodels[2],fid,"diagnostic","stokes");
73 _printf_("\n reading diagnostic hutter model data:\n");
74 CreateFemModel(&femmodels[3],fid,"diagnostic","hutter");
75 _printf_("\n reading surface and bed slope computation model data:\n");
76 CreateFemModel(&femmodels[4],fid,"slope_compute","");
77 _printf_("\n reading prognositc model data:\n");
78 CreateFemModel(&femmodels[5],fid,"prognostic","");
79
80 /*Do we run in 3d?, in which case we need thermal and melting also:*/
81 femmodels[0].parameters->FindParam((void*)&dim,"dim");
82 if(dim==3){
83 _printf_("read and create thermal finite element model:\n");
84 CreateFemModel(&femmodels[6],fid,"thermal","transient");
85 _printf_("read and create melting finite element model:\n");
86 CreateFemModel(&femmodels[7],fid,"melting","transient");
87 }
88
89 _printf_("initialize inputs:\n");
90 inputs=new ParameterInputs;
91 femmodels[5].parameters->FindParam((void*)&numberofnodes,"numberofnodes");
92
93 femmodels[5].parameters->FindParam((void*)&u_g,"u_g");
94 inputs->Add("velocity",u_g,3,numberofnodes);
95
96 femmodels[5].parameters->FindParam((void*)&m_g,"m_g");
97 inputs->Add("melting",m_g,1,numberofnodes);
98
99 femmodels[5].parameters->FindParam((void*)&a_g,"a_g");
100 inputs->Add("accumulation",a_g,1,numberofnodes);
101
102 femmodels[5].parameters->FindParam((void*)&dt,"dt");
103 inputs->Add("dt",dt);
104
105 _printf_("initialize results:\n");
106 results=new DataSet(ResultsEnum());
107
108 /*are we running the solution sequence, or a qmu wrapper around it? : */
109 femmodels[5].parameters->FindParam((void*)&qmu_analysis,"qmu_analysis");
110 if(!qmu_analysis){
111
112 /*run diagnostic analysis: */
113 _printf_("call computational core:\n");
114 transient_core(results,femmodels,inputs);
115 }
116 else{
117 /*run qmu analysis: */
118 _printf_("calling qmu analysis on transient core:\n");
119
120 #ifdef _HAVE_DAKOTA_
121 Qmux(&femmodels[0],inputs,TransientAnalysisEnum(),NoneAnalysisEnum());
122 #else
123 throw ErrorException(__FUNCT__," Dakota not present, cannot do qmu!");
124 #endif
125 }
126
127 /*Add analysis_type to results: */
128 result=new Result(results->Size()+1,0,1,"analysis_type","transient");
129 results->AddObject(result);
130
131 _printf_("process results:\n");
132 ProcessResults(&results,&femmodels[0],TransientAnalysisEnum());
133
134 _printf_("write results to disk:\n");
135 OutputResults(results,outputfilename);
136
137 _printf_("write lock file:\n");
138 femmodels[0].parameters->FindParam((void*)&waitonlock,"waitonlock");
139 if (waitonlock){
140 WriteLockFile(lockname);
141 }
142
143 _printf_("closing MPI and Petsc\n");
144 PetscFinalize();
145
146 /*end module: */
147 MODULEEND();
148
149 return 0; //unix success return;
150}
Note: See TracBrowser for help on using the repository browser.