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

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

New transient 2d solution in parallel

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