source: issm/trunk/src/c/solutions/diagnostic.cpp@ 4295

Last change on this file since 4295 was 4295, checked in by Mathieu Morlighem, 15 years ago

solution_type is now a SolutionEnum (not an AnalysisEnum) and some variables are now const for a better initialization

File size: 4.0 KB
RevLine 
[1]1/*!\file: diagnostic.cpp
2 * \brief: diagnostic 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"
[4236]13#include "../Container/Container.h"
[3751]14#include "../EnumDefinitions/EnumDefinitions.h"
[3775]15#include "../include/include.h"
[3913]16#include "../modules/modules.h"
[3895]17#include "./solutions.h"
[643]18
[1]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;
[3767]26 bool qmu_analysis=false;
27 bool control_analysis=false;
[4037]28 bool waitonlock=false;
[1]29
[3982]30 /*FemModel: */
31 FemModel* femmodel=NULL;
[643]32
[2489]33 /*time*/
34 double start, finish;
35 double start_core, finish_core;
36 double start_init, finish_init;
37
[4295]38 const int numanalyses=6;
39 const int analyses[numanalyses]={DiagnosticHorizAnalysisEnum,DiagnosticVertAnalysisEnum,DiagnosticStokesAnalysisEnum,DiagnosticHutterAnalysisEnum,SurfaceSlopeAnalysisEnum,BedSlopeAnalysisEnum};
40 const int solution_type=DiagnosticAnalysisEnum;
[4004]41
[1]42 MODULEBOOT();
43
44 #if !defined(_PARALLEL_) || (defined(_PARALLEL_) && !defined(_HAVE_PETSC_))
[3332]45 ISSMERROR(" parallel executable was compiled without support of parallel libraries!");
[1]46 #endif
47
[2489]48 /*Initialize Petsc and get start time*/
[1]49 PetscInitialize(&argc,&argv,(char *)0,"");
[2489]50 MPI_Barrier(MPI_COMM_WORLD); start=MPI_Wtime();
[1]51
52 /*Size and rank: */
53 MPI_Comm_rank(MPI_COMM_WORLD,&my_rank);
54 MPI_Comm_size(MPI_COMM_WORLD,&num_procs);
55
56 _printf_("recover , input file name and output file name:\n");
57 inputfilename=argv[2];
58 outputfilename=argv[3];
59 lockname=argv[4];
60
[2489]61 MPI_Barrier(MPI_COMM_WORLD); start_init=MPI_Wtime();
[1826]62
[1]63 /*Open handle to data on disk: */
[472]64 fid=pfopen(inputfilename,"rb");
[358]65
[4037]66 _printf_("create finite element model:\n");
[4295]67 femmodel=new FemModel(fid,solution_type,analyses,numanalyses);
[4055]68
69 /*add outputfilename in parameters: */
70 femmodel->parameters->AddObject(new StringParam(OutputFileNameEnum,outputfilename));
[1881]71
[1887]72 /*get parameters: */
[3982]73 femmodel->parameters->FindParam(&qmu_analysis,QmuAnalysisEnum);
74 femmodel->parameters->FindParam(&control_analysis,ControlAnalysisEnum);
75 femmodel->parameters->FindParam(&waitonlock,WaitOnLockEnum);
[1]76
[2489]77 MPI_Barrier(MPI_COMM_WORLD); finish_init=MPI_Wtime();
78
[732]79 /*are we running the solution sequence, or a qmu wrapper around it? : */
[586]80 if(!qmu_analysis){
[1911]81 if(!control_analysis){
[3887]82
[1911]83 _printf_("call computational core:\n");
[2489]84 MPI_Barrier(MPI_COMM_WORLD); start_core=MPI_Wtime( );
[4014]85 diagnostic_core(femmodel);
[2489]86 MPI_Barrier(MPI_COMM_WORLD); finish_core=MPI_Wtime( );
[2234]87
[2048]88 }
89 else{
[2234]90 /*run control analysis: */
91 _printf_("call computational core:\n");
[2489]92 MPI_Barrier(MPI_COMM_WORLD); start_core=MPI_Wtime( );
[4014]93 control_core(femmodel);
[2489]94 MPI_Barrier(MPI_COMM_WORLD); finish_core=MPI_Wtime( );
[2234]95
[2048]96 }
97
98 _printf_("write results to disk:\n");
[4211]99 OutputResultsx(femmodel->elements, femmodel->nodes, femmodel->vertices, femmodel->loads, femmodel->materials, femmodel->parameters,femmodel->results);
[586]100 }
101 else{
102 /*run qmu analysis: */
103 _printf_("calling qmu analysis on diagnostic core:\n");
[804]104
[662]105 #ifdef _HAVE_DAKOTA_
[2489]106 MPI_Barrier(MPI_COMM_WORLD); start_core=MPI_Wtime( );
[3982]107 Qmux(femmodel,DiagnosticAnalysisEnum,NoneAnalysisEnum);
[2489]108 MPI_Barrier(MPI_COMM_WORLD); finish_core=MPI_Wtime( );
[662]109 #else
[3332]110 ISSMERROR(" Dakota not present, cannot do qmu!");
[662]111 #endif
[586]112 }
[1]113
[2629]114 if (waitonlock>0){
115 _printf_("write lock file:\n");
[1]116 WriteLockFile(lockname);
117 }
[2048]118
[1911]119 /*Free ressources */
[3982]120 delete femmodel;
[1911]121
[2489]122 /*Get finish time and close*/
123 MPI_Barrier(MPI_COMM_WORLD); finish = MPI_Wtime( );
[3982]124 _printf_("\n %-34s %f seconds \n","FemModel initialization elapsed time:",finish_init-start_init);
[2489]125 _printf_(" %-34s %f seconds \n","Core solution elapsed time:",finish_core-start_core);
[3097]126 _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);
[2397]127 _printf_("closing MPI and Petsc\n");
128 PetscFinalize();
129
130 /*end module: */
131 MODULEEND();
132
[1]133 return 0; //unix success return;
134}
Note: See TracBrowser for help on using the repository browser.