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

Last change on this file since 1719 was 1719, checked in by seroussi, 16 years ago

changed cielo parallel timestep in yrs

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