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

Last change on this file since 3859 was 3859, checked in by seroussi, 15 years ago

some fixing for Dakota (still not working)

File size: 4.7 KB
RevLine 
[823]1/*!\file: transient.cpp
2 * \brief: transient 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
[3751]11#include "../objects/objects.h"
12#include "../shared/shared.h"
13#include "../DataSet/DataSet.h"
14#include "../EnumDefinitions/EnumDefinitions.h"
[3775]15#include "../include/include.h"
[3859]16#include "../modules.h"
[3751]17#include "./parallel.h"
[823]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 int numberofnodes;
[3767]27 bool qmu_analysis=false;
28 bool waitonlock=false;
[823]29
[1820]30 /*Model: */
31 Model* model=NULL;
[823]32 int dim=-1;
33
34 /*Results: */
35 DataSet* results=NULL;
[2112]36 DataSet* processed_results=NULL;
[1271]37 Result* result=NULL;
[823]38
39 Param* param=NULL;
40
[2492]41 /*time*/
42 double start, finish;
43 double start_core, finish_core;
44 double start_init, finish_init;
45
[823]46 MODULEBOOT();
47
48 #if !defined(_PARALLEL_) || (defined(_PARALLEL_) && !defined(_HAVE_PETSC_))
[3332]49 ISSMERROR(" parallel executable was compiled without support of parallel libraries!");
[823]50 #endif
51
[2492]52 /*Initialize Petsc and get start time*/
[823]53 PetscInitialize(&argc,&argv,(char *)0,"");
[2492]54 MPI_Barrier(MPI_COMM_WORLD); start=MPI_Wtime();
[823]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
65 /*Open handle to data on disk: */
66 fid=pfopen(inputfilename,"rb");
67
[1820]68 /*Initialize model structure: */
[2492]69 MPI_Barrier(MPI_COMM_WORLD); start_init=MPI_Wtime();
[1820]70 model=new Model();
71
[823]72 _printf_("read and create finite element model:\n");
73 _printf_("\n reading diagnostic horiz model data:\n");
[3567]74 model->AddFormulation(fid,DiagnosticAnalysisEnum,HorizAnalysisEnum);
[1881]75
[823]76 _printf_("\n reading diagnostic vert model data:\n");
[3567]77 model->AddFormulation(fid,DiagnosticAnalysisEnum,VertAnalysisEnum);
[1881]78
[823]79 _printf_("\n reading diagnostic stokes model data:\n");
[3567]80 model->AddFormulation(fid,DiagnosticAnalysisEnum,StokesAnalysisEnum);
[1881]81
[823]82 _printf_("\n reading diagnostic hutter model data:\n");
[3567]83 model->AddFormulation(fid,DiagnosticAnalysisEnum,HutterAnalysisEnum);
[1881]84
[823]85 _printf_("\n reading surface and bed slope computation model data:\n");
[3567]86 model->AddFormulation(fid,SlopecomputeAnalysisEnum);
[1881]87
[823]88 _printf_("\n reading prognositc model data:\n");
[3567]89 model->AddFormulation(fid,PrognosticAnalysisEnum);
[823]90
91 /*Do we run in 3d?, in which case we need thermal and melting also:*/
[3709]92 model->FindParam(&dim,DimEnum);
[823]93 if(dim==3){
94 _printf_("read and create thermal finite element model:\n");
[3567]95 model->AddFormulation(fid,ThermalAnalysisEnum,TransientAnalysisEnum);
[823]96 _printf_("read and create melting finite element model:\n");
[3567]97 model->AddFormulation(fid,MeltingAnalysisEnum,TransientAnalysisEnum);
[823]98 }
99
[1887]100 /*recover parameters: */
[3709]101 model->FindParam(&waitonlock,WaitOnLockEnum);
102 model->FindParam(&qmu_analysis,QmuAnalysisEnum);
[1887]103
[823]104 _printf_("initialize results:\n");
[3567]105 results=new DataSet(ResultsEnum);
[2492]106 MPI_Barrier(MPI_COMM_WORLD); finish_init=MPI_Wtime();
[823]107
108 /*are we running the solution sequence, or a qmu wrapper around it? : */
109 if(!qmu_analysis){
110
111 /*run diagnostic analysis: */
112 _printf_("call computational core:\n");
[2492]113 MPI_Barrier(MPI_COMM_WORLD); start_core=MPI_Wtime( );
[3751]114 transient_core(results,model);
[2492]115 MPI_Barrier(MPI_COMM_WORLD); finish_core=MPI_Wtime( );
[2048]116
117 /*Add analysis_type to results: */
118 result=new Result(results->Size()+1,0,1,"analysis_type","transient");
119 results->AddObject(result);
120
121 _printf_("process results:\n");
[3567]122 ProcessResults(&processed_results,results,model,TransientAnalysisEnum);
[2048]123
124 _printf_("write results to disk:\n");
[2113]125 OutputResults(processed_results,outputfilename);
[823]126 }
127 else{
128 /*run qmu analysis: */
129 _printf_("calling qmu analysis on transient core:\n");
130
131 #ifdef _HAVE_DAKOTA_
[2492]132 MPI_Barrier(MPI_COMM_WORLD); start_core=MPI_Wtime( );
[3751]133 Qmux(model,TransientAnalysisEnum,NoneAnalysisEnum);
[2492]134 MPI_Barrier(MPI_COMM_WORLD); finish_core=MPI_Wtime( );
[823]135 #else
[3332]136 ISSMERROR(" Dakota not present, cannot do qmu!");
[823]137 #endif
138 }
[2629]139
140 if (waitonlock>0){
141 _printf_("write lock file:\n");
[823]142 WriteLockFile(lockname);
143 }
[2397]144
[2013]145 /*Free ressources:*/
146 delete results;
[2112]147 delete processed_results;
[2013]148 delete model;
[2492]149
150 /*Get finish time and close*/
151 MPI_Barrier(MPI_COMM_WORLD); finish = MPI_Wtime( );
152 _printf_("\n %-34s %f seconds \n","Model initialization elapsed time:",finish_init-start_init);
153 _printf_(" %-34s %f seconds \n","Core solution elapsed time:",finish_core-start_core);
[3097]154 _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);
[823]155 _printf_("closing MPI and Petsc\n");
156 PetscFinalize();
157
158 /*end module: */
159 MODULEEND();
160
161 return 0; //unix success return;
162}
Note: See TracBrowser for help on using the repository browser.