xref: /phasta/phSolver/common/input_fform.cc (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen #include <fstream>
2*59599516SKenneth E. Jansen #include <stdio.h>
3*59599516SKenneth E. Jansen #include <stdlib.h>
4*59599516SKenneth E. Jansen #include <vector>
5*59599516SKenneth E. Jansen #include <string>
6*59599516SKenneth E. Jansen //MR CHANGE
7*59599516SKenneth E. Jansen #include <cstring>
8*59599516SKenneth E. Jansen //MR CHANGE END
9*59599516SKenneth E. Jansen 
10*59599516SKenneth E. Jansen #include "Input.h"
11*59599516SKenneth E. Jansen #include "common_c.h"
12*59599516SKenneth E. Jansen 
13*59599516SKenneth E. Jansen using namespace std; //::cout;
14*59599516SKenneth E. Jansen void print_error_code(int ierr);
15*59599516SKenneth E. Jansen int input_fform(char inpfname[]);
16*59599516SKenneth E. Jansen int SONFATH=0;
17*59599516SKenneth E. Jansen extern "C" char phasta_iotype[80];
18*59599516SKenneth E. Jansen 
19*59599516SKenneth E. Jansen extern "C" {
20*59599516SKenneth E. Jansen 
21*59599516SKenneth E. Jansen   int input_fform_(char inpfname[]) {
22*59599516SKenneth E. Jansen     return input_fform(inpfname);
23*59599516SKenneth E. Jansen   }
24*59599516SKenneth E. Jansen 
25*59599516SKenneth E. Jansen }
26*59599516SKenneth E. Jansen 
27*59599516SKenneth E. Jansen int input_fform(char inpfname[])
28*59599516SKenneth E. Jansen {
29*59599516SKenneth E. Jansen 
30*59599516SKenneth E. Jansen   int ierr = 0 ;
31*59599516SKenneth E. Jansen   int i,j;
32*59599516SKenneth E. Jansen   char* path_to_config = 0 ;
33*59599516SKenneth E. Jansen   char complete_filename[256];
34*59599516SKenneth E. Jansen 
35*59599516SKenneth E. Jansen   try {
36*59599516SKenneth E. Jansen     // get the input file stream
37*59599516SKenneth E. Jansen     path_to_config = getenv("PHASTA_CONFIG");
38*59599516SKenneth E. Jansen     if(path_to_config) strcpy(complete_filename, path_to_config);
39*59599516SKenneth E. Jansen     else strcpy(complete_filename,".");
40*59599516SKenneth E. Jansen     strcat(complete_filename, "/input.config");
41*59599516SKenneth E. Jansen     if(workfc.myrank==workfc.master) {
42*59599516SKenneth E. Jansen       printf("\n Complete Filename: %s \n", complete_filename);
43*59599516SKenneth E. Jansen       printf("\n Local Config: %s \n\n", inpfname);
44*59599516SKenneth E. Jansen     }
45*59599516SKenneth E. Jansen     string def(complete_filename);
46*59599516SKenneth E. Jansen     Input inp(inpfname,def);
47*59599516SKenneth E. Jansen 
48*59599516SKenneth E. Jansen #ifdef AMG
49*59599516SKenneth E. Jansen 
50*59599516SKenneth E. Jansen 	// AMG PARAMETERS
51*59599516SKenneth E. Jansen 
52*59599516SKenneth E. Jansen 	if ((string)inp.GetValue("Employ AMG") == "True" ) {
53*59599516SKenneth E. Jansen 
54*59599516SKenneth E. Jansen 	    amgvari.irun_amg = 1;
55*59599516SKenneth E. Jansen 
56*59599516SKenneth E. Jansen         amgvari.irun_amg_prec = inp.GetValue("Run AMG As CG-preconditioner");
57*59599516SKenneth E. Jansen 
58*59599516SKenneth E. Jansen         amgvarr.strong_eps       = inp.GetValue("Strong Criterion Eps");
59*59599516SKenneth E. Jansen 
60*59599516SKenneth E. Jansen         amgvarr.ramg_eps         = inp.GetValue("AMG Convergence Eps");
61*59599516SKenneth E. Jansen 
62*59599516SKenneth E. Jansen         amgvarr.ramg_relax       = inp.GetValue("AMG Relaxation Omega");
63*59599516SKenneth E. Jansen         amgvarr.ramg_trunc       = inp.GetValue("AMG Truncation Set");
64*59599516SKenneth E. Jansen 
65*59599516SKenneth E. Jansen         amgvari.iamg_verb        = inp.GetValue("AMG Verbosity");
66*59599516SKenneth E. Jansen 
67*59599516SKenneth E. Jansen         amgvari.iamg_neg_sten    = inp.GetValue("AMG Neg_Sten");
68*59599516SKenneth E. Jansen 
69*59599516SKenneth E. Jansen         amgvari.iamg_nlevel      = inp.GetValue("AMG Nlevel");
70*59599516SKenneth E. Jansen 
71*59599516SKenneth E. Jansen         amgvari.iamg_c_solver  = inp.GetValue("AMG Coarsest Solver");
72*59599516SKenneth E. Jansen 
73*59599516SKenneth E. Jansen         amgvari.iamg_init = 0;
74*59599516SKenneth E. Jansen         amgvari.iamg_setup_frez = inp.GetValue("AMG Freeze Setup");
75*59599516SKenneth E. Jansen         if ((string)inp.GetValue("AMG Interpolation Type")=="Standard")
76*59599516SKenneth E. Jansen             amgvari.iamg_interp = 1;
77*59599516SKenneth E. Jansen         else
78*59599516SKenneth E. Jansen             amgvari.iamg_interp = 0;
79*59599516SKenneth E. Jansen         amgvari.maxnev = inp.GetValue("AMG GGB nev");
80*59599516SKenneth E. Jansen         amgvari.maxncv = inp.GetValue("AMG GGB ncv");
81*59599516SKenneth E. Jansen 
82*59599516SKenneth E. Jansen         if ((string)inp.GetValue("AMG Smoother Type")=="Gauss-Seidel")
83*59599516SKenneth E. Jansen             amgvari.iamg_smoother = 1;
84*59599516SKenneth E. Jansen         else if ((string)inp.GetValue("AMG Smoother Type")=="ChebyShev")
85*59599516SKenneth E. Jansen             amgvari.iamg_smoother = 2;
86*59599516SKenneth E. Jansen         else if ((string)inp.GetValue("AMG Smoother Type")=="MLS")
87*59599516SKenneth E. Jansen             amgvari.iamg_smoother = 3;
88*59599516SKenneth E. Jansen         amgvarr.ramg_chebyratio       = inp.GetValue("AMG Chebyshev Eigenvalue ratio");
89*59599516SKenneth E. Jansen 
90*59599516SKenneth E. Jansen         amgvari.mlsdeg = inp.GetValue("AMG MLS Degree");
91*59599516SKenneth E. Jansen         amgvari.iamg_reduce = inp.GetValue("AMG Run Reduced Serial");
92*59599516SKenneth E. Jansen 	}
93*59599516SKenneth E. Jansen #endif
94*59599516SKenneth E. Jansen 
95*59599516SKenneth E. Jansen /////////////////////////////chen Sep 25 2009  Flow Control Parameters ////////
96*59599516SKenneth E. Jansen // Take BC from IC at inlet
97*59599516SKenneth E. Jansen       ctrlvari.iI2Binlet = (int)inp.GetValue("Take BC from IC at Inlet");
98*59599516SKenneth E. Jansen       if(ctrlvari.iI2Binlet > 0 )
99*59599516SKenneth E. Jansen         {
100*59599516SKenneth E. Jansen          ctrlvar.inletVelX = inp.GetValue("Inlet Bulk x Velocity");
101*59599516SKenneth E. Jansen         }
102*59599516SKenneth E. Jansen // set uniform outlet pressure
103*59599516SKenneth E. Jansen       ctrlvari.isetOutPres = (int)inp.GetValue("Set Outlet Pressure");
104*59599516SKenneth E. Jansen       if(ctrlvari.isetOutPres > 0 )
105*59599516SKenneth E. Jansen         {
106*59599516SKenneth E. Jansen          ctrlvar.outPres1 = inp.GetValue("Uniform Outlet Pressure");
107*59599516SKenneth E. Jansen         }
108*59599516SKenneth E. Jansen // set initial condition
109*59599516SKenneth E. Jansen       ctrlvari.isetInitial=(int)inp.GetValue("Specify Initial Conditions");
110*59599516SKenneth E. Jansen       if(ctrlvari.isetInitial==1)
111*59599516SKenneth E. Jansen         {
112*59599516SKenneth E. Jansen          ctrlvar.xvel_ini=inp.GetValue("Initial X Velocity");
113*59599516SKenneth E. Jansen          ctrlvar.yvel_ini=inp.GetValue("Initial Y Velocity");
114*59599516SKenneth E. Jansen          ctrlvar.zvel_ini=inp.GetValue("Initial Z Velocity");
115*59599516SKenneth E. Jansen          ctrlvar.temp_ini=inp.GetValue("Initial Temp");
116*59599516SKenneth E. Jansen          ctrlvar.pres_ini=inp.GetValue("Initial Pressure");
117*59599516SKenneth E. Jansen          ctrlvar.evis_ini=inp.GetValue("Initial Scalar 1");
118*59599516SKenneth E. Jansen         }
119*59599516SKenneth E. Jansen 
120*59599516SKenneth E. Jansen ///////////////////////////////////////////////////////////////////////////////
121*59599516SKenneth E. Jansen 
122*59599516SKenneth E. Jansen 
123*59599516SKenneth E. Jansen     // Disabled Features
124*59599516SKenneth E. Jansen 
125*59599516SKenneth E. Jansen     conpar.iALE = inp.GetValue("iALE");
126*59599516SKenneth E. Jansen     conpar.icoord = inp.GetValue("icoord");
127*59599516SKenneth E. Jansen     conpar.irs = inp.GetValue("irs");
128*59599516SKenneth E. Jansen     conpar.iexec = inp.GetValue("iexec");
129*59599516SKenneth E. Jansen     timpar.ntseq = inp.GetValue("ntseq");
130*59599516SKenneth E. Jansen     solpar.imap = inp.GetValue("imap");
131*59599516SKenneth E. Jansen 
132*59599516SKenneth E. Jansen 
133*59599516SKenneth E. Jansen     // Solution Control Keywords
134*59599516SKenneth E. Jansen 
135*59599516SKenneth E. Jansen     if((string)inp.GetValue("Equation of State") == "Incompressible") matdat.matflg[0][0] =-1 ;
136*59599516SKenneth E. Jansen     if((string)inp.GetValue("Equation of State") == "Compressible") matdat.matflg[0][0] =0;
137*59599516SKenneth E. Jansen     inpdat.Delt[0] = inp.GetValue("Time Step Size");
138*59599516SKenneth E. Jansen     inpdat.nstep[0] = inp.GetValue("Number of Timesteps");
139*59599516SKenneth E. Jansen     if((string)inp.GetValue("Viscous Control")=="Viscous") conpar.navier=1 ; else conpar.navier=0;
140*59599516SKenneth E. Jansen 
141*59599516SKenneth E. Jansen     if ((string)inp.GetValue("Turbulence Model") == "No-Model" ) {
142*59599516SKenneth E. Jansen       turbvari.irans = 0;
143*59599516SKenneth E. Jansen       turbvari.iles  = 0;
144*59599516SKenneth E. Jansen     } else if ((string)inp.GetValue("Turbulence Model") == "LES" ) {
145*59599516SKenneth E. Jansen       turbvari.iles  = 1;
146*59599516SKenneth E. Jansen       turbvari.irans = 0;
147*59599516SKenneth E. Jansen     } else if ((string)inp.GetValue("Turbulence Model") == "RANS-SA" ) {
148*59599516SKenneth E. Jansen       turbvari.iles  = 0;
149*59599516SKenneth E. Jansen       turbvari.irans = -1;
150*59599516SKenneth E. Jansen     } else if ((string)inp.GetValue("Turbulence Model") == "RANS" ) {
151*59599516SKenneth E. Jansen       turbvari.iles  = 0;
152*59599516SKenneth E. Jansen       turbvari.irans = -1; // assume S-A for backward compatibility
153*59599516SKenneth E. Jansen     } else if ((string)inp.GetValue("Turbulence Model") == "RANS-KE" ) {
154*59599516SKenneth E. Jansen       turbvari.iles  = 0;
155*59599516SKenneth E. Jansen       turbvari.irans = -2;
156*59599516SKenneth E. Jansen     } else if ((string)inp.GetValue("Turbulence Model") == "DES97" ) {
157*59599516SKenneth E. Jansen       turbvari.iles  = -1;
158*59599516SKenneth E. Jansen       turbvari.irans = -1;
159*59599516SKenneth E. Jansen     } else if ((string)inp.GetValue("Turbulence Model") == "DDES" ) {
160*59599516SKenneth E. Jansen       turbvari.iles  = -2;
161*59599516SKenneth E. Jansen       turbvari.irans = -1;
162*59599516SKenneth E. Jansen 
163*59599516SKenneth E. Jansen     } else {
164*59599516SKenneth E. Jansen       cout << " Turbulence Model: Only Legal Values ( No-Model, LES, RANS-SA, RANS-KE, DES97, DDES )";
165*59599516SKenneth E. Jansen       cout << endl;
166*59599516SKenneth E. Jansen       exit(1);
167*59599516SKenneth E. Jansen     }
168*59599516SKenneth E. Jansen 
169*59599516SKenneth E. Jansen  //   if (turbvari.iles*turbvari.irans!=0) turbvar.eles=
170*59599516SKenneth E. Jansen  //                                          inp.GetValue("DES Edge Length");
171*59599516SKenneth E. Jansen 
172*59599516SKenneth E. Jansen     if (turbvari.irans<0 && turbvari.iles<0)
173*59599516SKenneth E. Jansen       turbvar.DES_SA_hmin=(double)inp.GetValue("DES SA Minimum Edge Length");
174*59599516SKenneth E. Jansen 
175*59599516SKenneth E. Jansen     int solflow, solheat , solscalr, ilset;
176*59599516SKenneth E. Jansen     ((string)inp.GetValue("Solve Flow") == "True")? solflow=1:solflow=0;
177*59599516SKenneth E. Jansen     ((string)inp.GetValue("Solve Heat") == "True")? solheat=1:solheat=0;
178*59599516SKenneth E. Jansen     //for compressible solheat= False so
179*59599516SKenneth E. Jansen     if((string)inp.GetValue("Equation of State") == "Compressible") solheat=0;
180*59599516SKenneth E. Jansen     ilset = (int)inp.GetValue("Solve Level Set");
181*59599516SKenneth E. Jansen     solscalr = (int)inp.GetValue("Solve Scalars");
182*59599516SKenneth E. Jansen     solscalr += ilset;
183*59599516SKenneth E. Jansen     if(turbvari.irans == -1) solscalr++;
184*59599516SKenneth E. Jansen     if(turbvari.irans == -2) solscalr=solscalr+2;
185*59599516SKenneth E. Jansen     if ( solscalr > 4 ) {
186*59599516SKenneth E. Jansen       cout << " Only Four Scalars are supported \n";
187*59599516SKenneth E. Jansen       cout <<" Please reduce number of scalars \n";
188*59599516SKenneth E. Jansen       exit(1);
189*59599516SKenneth E. Jansen     }
190*59599516SKenneth E. Jansen     inpdat.impl[0] = 10*solflow+solscalr*100+solheat;
191*59599516SKenneth E. Jansen 
192*59599516SKenneth E. Jansen     levlset.iLSet = ilset;
193*59599516SKenneth E. Jansen     if( ilset > 0) {
194*59599516SKenneth E. Jansen     levlset.epsilon_ls = inp.GetValue("Number of Elements Across Interface");
195*59599516SKenneth E. Jansen     levlset.epsilon_lsd = inp.GetValue("Number of Elements Across Interface for Redistancing");
196*59599516SKenneth E. Jansen     levlset.dtlset = inp.GetValue("Pseudo Time step for Redistancing");
197*59599516SKenneth E. Jansen     levlset.iExpLSSclr2 = inp.GetValue("Explicit Solve for Redistance Field");
198*59599516SKenneth E. Jansen     levlset.iExpLSSclr1 = inp.GetValue("Explicit Solve for Scalar 1 Field");
199*59599516SKenneth E. Jansen     if ((string)inp.GetValue("Apply Volume Constraint") == "True" ) {
200*59599516SKenneth E. Jansen       levlset.ivconstraint = 1; }
201*59599516SKenneth E. Jansen     else if((string)inp.GetValue("Apply Volume Constraint") == "False" ) {
202*59599516SKenneth E. Jansen       levlset.ivconstraint = 0; }
203*59599516SKenneth E. Jansen     else {
204*59599516SKenneth E. Jansen       cout << "Apply Volume Constraint: Only Legal Values (True, False) ";
205*59599516SKenneth E. Jansen       cout << endl;
206*59599516SKenneth E. Jansen       exit(1);
207*59599516SKenneth E. Jansen     }
208*59599516SKenneth E. Jansen     }
209*59599516SKenneth E. Jansen 
210*59599516SKenneth E. Jansen     vector<double> vec;
211*59599516SKenneth E. Jansen 
212*59599516SKenneth E. Jansen     // OUTPUT CONTROL KEY WORDS.
213*59599516SKenneth E. Jansen 
214*59599516SKenneth E. Jansen     conpar.necho = inp.GetValue("Verbosity Level");
215*59599516SKenneth E. Jansen     outpar.ntout = inp.GetValue("Number of Timesteps between Restarts");
216*59599516SKenneth E. Jansen     if((string)inp.GetValue("Print Statistics") == "True") outpar.ioform = 2;
217*59599516SKenneth E. Jansen     else outpar.ioform = 1;
218*59599516SKenneth E. Jansen 
219*59599516SKenneth E. Jansen     if((string)inp.GetValue("Print Wall Fluxes") == "True") outpar.iowflux = 1;
220*59599516SKenneth E. Jansen     else outpar.iowflux = 0;
221*59599516SKenneth E. Jansen 
222*59599516SKenneth E. Jansen     if((string)inp.GetValue("Print FieldView") == "True") outpar.iofieldv = 1;
223*59599516SKenneth E. Jansen     else outpar.iofieldv = 0;
224*59599516SKenneth E. Jansen 
225*59599516SKenneth E. Jansen     if((string)inp.GetValue("Print ybar") == "True") outpar.ioybar = 1;
226*59599516SKenneth E. Jansen     else outpar.ioybar = 0;
227*59599516SKenneth E. Jansen 
228*59599516SKenneth E. Jansen     if((string)inp.GetValue("Print vorticity") == "True") outpar.ivort = 1;
229*59599516SKenneth E. Jansen     else outpar.ivort = 0;
230*59599516SKenneth E. Jansen 
231*59599516SKenneth E. Jansen     outpar.nstepsincycle = inp.GetValue("Number of Steps in a Cycle");
232*59599516SKenneth E. Jansen     outpar.nphasesincycle = inp.GetValue("Number of Phases in a Cycle");
233*59599516SKenneth E. Jansen     outpar.ncycles_startphaseavg = inp.GetValue("Number of Initial Cycles to Skip in Phase Average");
234*59599516SKenneth E. Jansen 
235*59599516SKenneth E. Jansen     strcpy( outpar.iotype , ((string)inp.GetValue("Data Block Format")).c_str());
236*59599516SKenneth E. Jansen     strcpy( phasta_iotype , ((string)inp.GetValue("Data Block Format")).c_str());
237*59599516SKenneth E. Jansen     SONFATH = inp.GetValue("Number of Father Nodes");
238*59599516SKenneth E. Jansen 
239*59599516SKenneth E. Jansen     if((string)inp.GetValue("Print Residual at End of Step") == "True") genpar.lstres = 1;
240*59599516SKenneth E. Jansen     else genpar.lstres = 0;
241*59599516SKenneth E. Jansen 
242*59599516SKenneth E. Jansen     if((string)inp.GetValue("Print Error Indicators") == "True") turbvar.ierrcalc = 1;
243*59599516SKenneth E. Jansen     else turbvar.ierrcalc = 0;
244*59599516SKenneth E. Jansen 
245*59599516SKenneth E. Jansen     if((string)inp.GetValue("Print Velocity Hessian") == "True") turbvar.ihessian = 1;
246*59599516SKenneth E. Jansen     else turbvar.ihessian = 0;
247*59599516SKenneth E. Jansen 
248*59599516SKenneth E. Jansen     if ( turbvar.ierrcalc == 1 )
249*59599516SKenneth E. Jansen         turbvari.ierrsmooth = inp.GetValue("Number of Error Smoothing Iterations");
250*59599516SKenneth E. Jansen 
251*59599516SKenneth E. Jansen     for(i=0;i<MAXSURF+1; i++) aerfrc.nsrflist[i] = 0;
252*59599516SKenneth E. Jansen     int nsrfCM = inp.GetValue("Number of Force Surfaces");
253*59599516SKenneth E. Jansen     if (nsrfCM > 0) {
254*59599516SKenneth E. Jansen       vector<int> ivec = inp.GetValue("Surface ID's for Force Calculation");
255*59599516SKenneth E. Jansen       for(i=0; i< nsrfCM; i++){
256*59599516SKenneth E. Jansen         aerfrc.nsrflist[ivec[i]] = 1;
257*59599516SKenneth E. Jansen         //        cout <<"surface in force list "<< ivec[i] << endl;
258*59599516SKenneth E. Jansen       }
259*59599516SKenneth E. Jansen       ivec.erase(ivec.begin(),ivec.end());
260*59599516SKenneth E. Jansen     }
261*59599516SKenneth E. Jansen 
262*59599516SKenneth E. Jansen     aerfrc.isrfIM = inp.GetValue("Surface ID for Integrated Mass");
263*59599516SKenneth E. Jansen 
264*59599516SKenneth E. Jansen     outpar.iv_rankpercore = inp.GetValue("Ranks per core");
265*59599516SKenneth E. Jansen     outpar.iv_corepernode = inp.GetValue("Cores per node");
266*59599516SKenneth E. Jansen 
267*59599516SKenneth E. Jansen     turbvari.iramp=0;
268*59599516SKenneth E. Jansen     if((string)inp.GetValue("Ramp Inflow") == "True") turbvari.iramp=1;
269*59599516SKenneth E. Jansen     if(turbvari.iramp == 1) {
270*59599516SKenneth E. Jansen 	vec = inp.GetValue("Mdot Ramp Inflow Start and Stop");
271*59599516SKenneth E. Jansen     	for(i=0; i<2 ; i++){
272*59599516SKenneth E. Jansen         	turbvar.rampmdot[0][i]=vec[i];
273*59599516SKenneth E. Jansen     	}
274*59599516SKenneth E. Jansen     	vec = inp.GetValue("Mdot Ramp Lower FC Start and Stop");
275*59599516SKenneth E. Jansen     	for(i=0; i<2 ; i++){
276*59599516SKenneth E. Jansen          	turbvar.rampmdot[1][i]=vec[i];
277*59599516SKenneth E. Jansen     	}
278*59599516SKenneth E. Jansen     	vec = inp.GetValue("Mdot Ramp Upper FC Start and Stop");
279*59599516SKenneth E. Jansen     	for(i=0; i<2 ; i++){
280*59599516SKenneth E. Jansen         	turbvar.rampmdot[2][i]=vec[i];
281*59599516SKenneth E. Jansen     	}
282*59599516SKenneth E. Jansen     }
283*59599516SKenneth E. Jansen 
284*59599516SKenneth E. Jansen //Limiting
285*59599516SKenneth E. Jansen     vec = inp.GetValue("Limit u1");
286*59599516SKenneth E. Jansen     for(i=0; i<3 ; i++){
287*59599516SKenneth E. Jansen       turbvar.ylimit[0][i] = vec[i];
288*59599516SKenneth E. Jansen     }
289*59599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
290*59599516SKenneth E. Jansen 
291*59599516SKenneth E. Jansen     vec = inp.GetValue("Limit u2");
292*59599516SKenneth E. Jansen     for(i=0; i<3 ; i++){
293*59599516SKenneth E. Jansen       turbvar.ylimit[1][i] = vec[i];
294*59599516SKenneth E. Jansen     }
295*59599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
296*59599516SKenneth E. Jansen 
297*59599516SKenneth E. Jansen     vec = inp.GetValue("Limit u3");
298*59599516SKenneth E. Jansen     for(i=0; i<3 ; i++){
299*59599516SKenneth E. Jansen       turbvar.ylimit[2][i] = vec[i];
300*59599516SKenneth E. Jansen     }
301*59599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
302*59599516SKenneth E. Jansen 
303*59599516SKenneth E. Jansen     vec = inp.GetValue("Limit Pressure");
304*59599516SKenneth E. Jansen     for(i=0; i<3 ; i++){
305*59599516SKenneth E. Jansen       turbvar.ylimit[3][i] = vec[i];
306*59599516SKenneth E. Jansen     }
307*59599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
308*59599516SKenneth E. Jansen 
309*59599516SKenneth E. Jansen     vec = inp.GetValue("Limit Temperature");
310*59599516SKenneth E. Jansen     for(i=0; i<3 ; i++){
311*59599516SKenneth E. Jansen       turbvar.ylimit[4][i] = vec[i];
312*59599516SKenneth E. Jansen     }
313*59599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
314*59599516SKenneth E. Jansen 
315*59599516SKenneth E. Jansen     //Material Properties Keywords
316*59599516SKenneth E. Jansen     matdat.nummat = levlset.iLSet+1;
317*59599516SKenneth E. Jansen     if((string)inp.GetValue("Shear Law") == "Constant Viscosity")
318*59599516SKenneth E. Jansen       for(i=0; i < levlset.iLSet+1; i++) matdat.matflg[i][1] = 0;
319*59599516SKenneth E. Jansen 
320*59599516SKenneth E. Jansen     if((string)inp.GetValue("Bulk Viscosity Law") == "Constant Bulk Viscosity")
321*59599516SKenneth E. Jansen       for(i=0; i < levlset.iLSet+1; i++) matdat.matflg[i][2] = 0;
322*59599516SKenneth E. Jansen 
323*59599516SKenneth E. Jansen     mmatpar.pr = inp.GetValue("Prandtl Number");
324*59599516SKenneth E. Jansen 
325*59599516SKenneth E. Jansen     if((string)inp.GetValue("Conductivity Law") == "Constant Conductivity")
326*59599516SKenneth E. Jansen       for(i=0; i < levlset.iLSet+1; i++) matdat.matflg[i][3] = 0;
327*59599516SKenneth E. Jansen 
328*59599516SKenneth E. Jansen     vec = inp.GetValue("Density");
329*59599516SKenneth E. Jansen     for(i=0; i< levlset.iLSet +1 ; i++){
330*59599516SKenneth E. Jansen       matdat.datmat[i][0][0] = vec[i];
331*59599516SKenneth E. Jansen     }
332*59599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
333*59599516SKenneth E. Jansen 
334*59599516SKenneth E. Jansen     vec = inp.GetValue("Viscosity");
335*59599516SKenneth E. Jansen     for(i=0; i< levlset.iLSet +1 ; i++){
336*59599516SKenneth E. Jansen       matdat.datmat[i][1][0] = vec[i];
337*59599516SKenneth E. Jansen     }
338*59599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
339*59599516SKenneth E. Jansen 
340*59599516SKenneth E. Jansen //      vec = inp.GetValue("Specific Heat");
341*59599516SKenneth E. Jansen     for(i=0; i< levlset.iLSet +1 ; i++){
342*59599516SKenneth E. Jansen       matdat.datmat[i][2][0] = 0;
343*59599516SKenneth E. Jansen     }
344*59599516SKenneth E. Jansen //      vec.erase(vec.begin(),vec.end());
345*59599516SKenneth E. Jansen 
346*59599516SKenneth E. Jansen     vec = inp.GetValue("Thermal Conductivity");
347*59599516SKenneth E. Jansen     for(i=0; i< levlset.iLSet +1 ; i++){
348*59599516SKenneth E. Jansen       matdat.datmat[i][3][0] = vec[i];
349*59599516SKenneth E. Jansen     }
350*59599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
351*59599516SKenneth E. Jansen 
352*59599516SKenneth E. Jansen     vec = inp.GetValue("Scalar Diffusivity");
353*59599516SKenneth E. Jansen     for(i=0; i< solscalr ; i++){
354*59599516SKenneth E. Jansen       sclrs.scdiff[i] = vec[i];
355*59599516SKenneth E. Jansen     }
356*59599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
357*59599516SKenneth E. Jansen 
358*59599516SKenneth E. Jansen     if((string)inp.GetValue("Zero Mean Pressure") == "True")
359*59599516SKenneth E. Jansen       turbvar.pzero=1;
360*59599516SKenneth E. Jansen 
361*59599516SKenneth E. Jansen     turbvar.rmutarget = inp.GetValue("Target Viscosity For Step NSTEP");
362*59599516SKenneth E. Jansen 
363*59599516SKenneth E. Jansen     if ( (string)inp.GetValue("Body Force Option") == "None" ) {
364*59599516SKenneth E. Jansen       for( i=0; i< levlset.iLSet +1 ; i++)  matdat.matflg[i][4] = 0;
365*59599516SKenneth E. Jansen     }
366*59599516SKenneth E. Jansen     else if ( (string)inp.GetValue("Body Force Option") == "Vector" ) {
367*59599516SKenneth E. Jansen       for( i=0; i< levlset.iLSet +1 ; i++)  matdat.matflg[i][4] = 1;
368*59599516SKenneth E. Jansen     }
369*59599516SKenneth E. Jansen     else if ( (string)inp.GetValue("Body Force Option") == "User e3source.f" ) {
370*59599516SKenneth E. Jansen       for( i=0; i< levlset.iLSet +1 ; i++) matdat.matflg[i][4] = 3;
371*59599516SKenneth E. Jansen     }
372*59599516SKenneth E. Jansen     else if ( (string)inp.GetValue("Body Force Option") == "Boussinesq" ) {
373*59599516SKenneth E. Jansen       for(i=0; i< levlset.iLSet +1 ; i++) matdat.matflg[i][4] = 2;
374*59599516SKenneth E. Jansen     }
375*59599516SKenneth E. Jansen     else if ( (string)inp.GetValue("Body Force Option") == "Cooling Analytic" ) {
376*59599516SKenneth E. Jansen       for(i=0; i< levlset.iLSet +1 ; i++) matdat.matflg[i][4] = 4;
377*59599516SKenneth E. Jansen     }
378*59599516SKenneth E. Jansen     else if ( (string)inp.GetValue("Body Force Option") == "Cooling Initial Condition" ) {
379*59599516SKenneth E. Jansen       for(i=0; i< levlset.iLSet +1 ; i++) matdat.matflg[i][4] = 5;
380*59599516SKenneth E. Jansen     }
381*59599516SKenneth E. Jansen 
382*59599516SKenneth E. Jansen     // the following block of stuff is common to all cooling type sponges.
383*59599516SKenneth E. Jansen     // Specific stuff belongs in the conditionals above
384*59599516SKenneth E. Jansen 
385*59599516SKenneth E. Jansen     if(matdat.matflg[0][4] >=4) {
386*59599516SKenneth E. Jansen       spongevar.betamax = inp.GetValue("Maximum Value of Sponge Parameter");
387*59599516SKenneth E. Jansen       spongevar.zinsponge = inp.GetValue("Inflow Cooling Sponge Ends at z");
388*59599516SKenneth E. Jansen       spongevar.zoutsponge= inp.GetValue("Outflow Cooling Sponge Begins at z");
389*59599516SKenneth E. Jansen       spongevar.radsponge = inp.GetValue("Radial Cooling Sponge Begins at r");
390*59599516SKenneth E. Jansen       spongevar.grthosponge = inp.GetValue("Sponge Growth Coefficient Outflow");
391*59599516SKenneth E. Jansen       spongevar.grthisponge = inp.GetValue("Sponge Growth Coefficient Inflow");
392*59599516SKenneth E. Jansen 
393*59599516SKenneth E. Jansen 
394*59599516SKenneth E. Jansen       spongevar.spongecontinuity = 0;
395*59599516SKenneth E. Jansen       spongevar.spongemomentum1 = 0;
396*59599516SKenneth E. Jansen       spongevar.spongemomentum2 = 0;
397*59599516SKenneth E. Jansen       spongevar.spongemomentum3 = 0;
398*59599516SKenneth E. Jansen       spongevar.spongeenergy = 0;
399*59599516SKenneth E. Jansen 
400*59599516SKenneth E. Jansen       if((string)inp.GetValue("Sponge for Continuity Equation") == "True")
401*59599516SKenneth E. Jansen 	spongevar.spongecontinuity = 1;
402*59599516SKenneth E. Jansen       if((string)inp.GetValue("Sponge for x Momentum Equation") == "True")
403*59599516SKenneth E. Jansen 	spongevar.spongemomentum1 = 1;
404*59599516SKenneth E. Jansen       if((string)inp.GetValue("Sponge for y Momentum Equation") == "True")
405*59599516SKenneth E. Jansen 	spongevar.spongemomentum2 = 1;
406*59599516SKenneth E. Jansen       if((string)inp.GetValue("Sponge for z Momentum Equation") == "True")
407*59599516SKenneth E. Jansen 	spongevar.spongemomentum3 = 1;
408*59599516SKenneth E. Jansen       if((string)inp.GetValue("Sponge for Energy Equation") == "True")
409*59599516SKenneth E. Jansen 	spongevar.spongeenergy = 1;
410*59599516SKenneth E. Jansen 
411*59599516SKenneth E. Jansen     }
412*59599516SKenneth E. Jansen 
413*59599516SKenneth E. Jansen     vec = inp.GetValue("Body Force");
414*59599516SKenneth E. Jansen     for(i=0; i< levlset.iLSet +1 ; i++){
415*59599516SKenneth E. Jansen       matdat.datmat[i][4][0] = vec[0+i*3];
416*59599516SKenneth E. Jansen       matdat.datmat[i][4][1] = vec[1+i*3];
417*59599516SKenneth E. Jansen       matdat.datmat[i][4][2] = vec[2+i*3];
418*59599516SKenneth E. Jansen     }
419*59599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
420*59599516SKenneth E. Jansen 
421*59599516SKenneth E. Jansen     vec = inp.GetValue("Body Force Pressure Gradient");
422*59599516SKenneth E. Jansen     for(i=0; i< levlset.iLSet +1 ; i++){
423*59599516SKenneth E. Jansen       matdat.datmat[i][6][0] = vec[0+i*3];
424*59599516SKenneth E. Jansen       matdat.datmat[i][6][1] = vec[1+i*3];
425*59599516SKenneth E. Jansen       matdat.datmat[i][6][2] = vec[2+i*3];
426*59599516SKenneth E. Jansen     }
427*59599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
428*59599516SKenneth E. Jansen 
429*59599516SKenneth E. Jansen     if ( (string)inp.GetValue("Surface Tension Option") == "No" ){
430*59599516SKenneth E. Jansen         genpar.isurf = 0;
431*59599516SKenneth E. Jansen     }
432*59599516SKenneth E. Jansen     else if ((string)inp.GetValue("Surface Tension Option") == "Yes" ){
433*59599516SKenneth E. Jansen         genpar.isurf = 1;
434*59599516SKenneth E. Jansen     }
435*59599516SKenneth E. Jansen     else {
436*59599516SKenneth E. Jansen       cout << " Surface Tension: Only Legal Values (Yes, No) ";
437*59599516SKenneth E. Jansen       cout << endl;
438*59599516SKenneth E. Jansen       exit(1);
439*59599516SKenneth E. Jansen     }
440*59599516SKenneth E. Jansen     if( genpar.isurf > 0) {
441*59599516SKenneth E. Jansen       genpar.Bo = inp.GetValue("Bond Number");
442*59599516SKenneth E. Jansen     }
443*59599516SKenneth E. Jansen 
444*59599516SKenneth E. Jansen     genpar.EntropyPressure = inp.GetValue("Entropy Form of Pressure Constraint on Weight Space");
445*59599516SKenneth E. Jansen 
446*59599516SKenneth E. Jansen 
447*59599516SKenneth E. Jansen     if ( (string)inp.GetValue("Rotating Frame of Reference") == "True" ) {
448*59599516SKenneth E. Jansen       matdat.matflg[0][5] = 1;
449*59599516SKenneth E. Jansen       vec = inp.GetValue("Rotating Frame of Reference Rotation Rate");
450*59599516SKenneth E. Jansen       matdat.datmat[0][5][0] = vec[0];
451*59599516SKenneth E. Jansen       matdat.datmat[0][5][1] = vec[1];
452*59599516SKenneth E. Jansen       matdat.datmat[0][5][2] = vec[2];
453*59599516SKenneth E. Jansen       vec.erase(vec.begin(),vec.end());
454*59599516SKenneth E. Jansen     }
455*59599516SKenneth E. Jansen     else {
456*59599516SKenneth E. Jansen       matdat.matflg[0][5] = 0;
457*59599516SKenneth E. Jansen       matdat.datmat[0][5][0] = 0.;
458*59599516SKenneth E. Jansen       matdat.datmat[0][5][1] = 0.;
459*59599516SKenneth E. Jansen       matdat.datmat[0][5][2] = 0.;
460*59599516SKenneth E. Jansen     }
461*59599516SKenneth E. Jansen 
462*59599516SKenneth E. Jansen 
463*59599516SKenneth E. Jansen     //Linear Solver parameters
464*59599516SKenneth E. Jansen     if( (string)inp.GetValue("Solver Type") =="ACUSIM with P Projection" ){
465*59599516SKenneth E. Jansen       incomp.iprjFlag = 0; incomp.ipresPrjFlag=1;}
466*59599516SKenneth E. Jansen     else if ( (string)inp.GetValue("Solver Type") =="ACUSIM" ){
467*59599516SKenneth E. Jansen       incomp.iprjFlag = 0; incomp.ipresPrjFlag=0;}
468*59599516SKenneth E. Jansen     else if( (string)inp.GetValue("Solver Type") =="ACUSIM with Velocity Projection" ){
469*59599516SKenneth E. Jansen       incomp.iprjFlag = 1; incomp.ipresPrjFlag=0;}
470*59599516SKenneth E. Jansen     else if( (string)inp.GetValue("Solver Type") =="ACUSIM with Full Projection" ){
471*59599516SKenneth E. Jansen       incomp.iprjFlag = 1; incomp.ipresPrjFlag=1;}
472*59599516SKenneth E. Jansen     else if( (string)inp.GetValue("Solver Type") =="GMRES Matrix Free"){
473*59599516SKenneth E. Jansen       inpdat.impl[0] += 10*solflow;}
474*59599516SKenneth E. Jansen     else if( (string)inp.GetValue("Solver Type") =="GMRES EBE"){
475*59599516SKenneth E. Jansen       inpdat.impl[0] += 20*solflow;}
476*59599516SKenneth E. Jansen     //GMRES sparse is assumed default and has the value of 10, MFG 20,
477*59599516SKenneth E. Jansen     // EBE 30
478*59599516SKenneth E. Jansen 
479*59599516SKenneth E. Jansen 
480*59599516SKenneth E. Jansen     //    inpdat.niter[0] = inp.GetValue("Number of Solves per Time Step");
481*59599516SKenneth E. Jansen     solpar.nGMRES = inp.GetValue("Number of GMRES Sweeps per Solve");
482*59599516SKenneth E. Jansen     solpar.Kspace = inp.GetValue("Number of Krylov Vectors per GMRES Sweep");
483*59599516SKenneth E. Jansen     inpdat.LHSupd[0] = inp.GetValue("Number of Solves per Left-hand-side Formation");
484*59599516SKenneth E. Jansen     inpdat.epstol[0] = inp.GetValue("Tolerance on Momentum Equations");
485*59599516SKenneth E. Jansen     incomp.prestol = inp.GetValue("Tolerance on ACUSIM Pressure Projection");
486*59599516SKenneth E. Jansen     incomp.minIters = inp.GetValue("Minimum Number of Iterations per Nonlinear Iteration");
487*59599516SKenneth E. Jansen     incomp.maxIters = inp.GetValue("Maximum Number of Iterations per Nonlinear Iteration");
488*59599516SKenneth E. Jansen     inpdat.deltol[0][0]=inp.GetValue("Velocity Delta Ratio");
489*59599516SKenneth E. Jansen     inpdat.deltol[1][0]=inp.GetValue("Pressure Delta Ratio");
490*59599516SKenneth E. Jansen     incomp.nPrjs = inp.GetValue("Number of Velocity Projection Vectors");
491*59599516SKenneth E. Jansen     incomp.nPresPrjs = inp.GetValue("Number of Pressure Projection Vectors");
492*59599516SKenneth E. Jansen     incomp.iverbose = inp.GetValue("ACUSIM Verbosity Level");
493*59599516SKenneth E. Jansen 
494*59599516SKenneth E. Jansen     if(solheat==1){
495*59599516SKenneth E. Jansen       inpdat.epstol[1]=inp.GetValue("Temperature Solver Tolerance");
496*59599516SKenneth E. Jansen       inpdat.LHSupd[1]=inp.GetValue("Number of Solves of Temperature per Left-hand-side Formation");
497*59599516SKenneth E. Jansen     }
498*59599516SKenneth E. Jansen 
499*59599516SKenneth E. Jansen     // The following is where you should put any inputs that are able to
500*59599516SKenneth E. Jansen     // input differently for each scalar.  It is a little tedious in the code
501*59599516SKenneth E. Jansen     // but it should make the solver.inp easier to understand. Note this will
502*59599516SKenneth E. Jansen     // require some care with regression tests.
503*59599516SKenneth E. Jansen 
504*59599516SKenneth E. Jansen 
505*59599516SKenneth E. Jansen     if(solscalr>0){
506*59599516SKenneth E. Jansen       inpdat.epstol[2]=inp.GetValue("Scalar 1 Solver Tolerance");
507*59599516SKenneth E. Jansen       inpdat.LHSupd[2]=inp.GetValue("Number of Solves of Scalar 1 per Left-hand-side Formation");
508*59599516SKenneth E. Jansen 
509*59599516SKenneth E. Jansen       vec = inp.GetValue("Limit Scalar 1");
510*59599516SKenneth E. Jansen       for(i=0; i<3 ; i++){
511*59599516SKenneth E. Jansen         turbvar.ylimit[5][i] = vec[i];
512*59599516SKenneth E. Jansen       }
513*59599516SKenneth E. Jansen       vec.erase(vec.begin(),vec.end());
514*59599516SKenneth E. Jansen     }
515*59599516SKenneth E. Jansen 
516*59599516SKenneth E. Jansen     if(solscalr>1){
517*59599516SKenneth E. Jansen       inpdat.epstol[3]=inp.GetValue("Scalar 2 Solver Tolerance");
518*59599516SKenneth E. Jansen       inpdat.LHSupd[3]=inp.GetValue("Number of Solves of Scalar 2 per Left-hand-side Formation");
519*59599516SKenneth E. Jansen 
520*59599516SKenneth E. Jansen       vec = inp.GetValue("Limit Scalar 2");
521*59599516SKenneth E. Jansen       for(i=0; i<3 ; i++){
522*59599516SKenneth E. Jansen         turbvar.ylimit[6][i] = vec[i];
523*59599516SKenneth E. Jansen       }
524*59599516SKenneth E. Jansen       vec.erase(vec.begin(),vec.end());
525*59599516SKenneth E. Jansen     }
526*59599516SKenneth E. Jansen 
527*59599516SKenneth E. Jansen     if(solscalr>2){
528*59599516SKenneth E. Jansen       inpdat.epstol[4]=inp.GetValue("Scalar 3 Solver Tolerance");
529*59599516SKenneth E. Jansen       inpdat.LHSupd[4]=inp.GetValue("Number of Solves of Scalar 3 per Left-hand-side Formation");
530*59599516SKenneth E. Jansen 
531*59599516SKenneth E. Jansen       vec = inp.GetValue("Limit Scalar 3");
532*59599516SKenneth E. Jansen       for(i=0; i<3 ; i++){
533*59599516SKenneth E. Jansen         turbvar.ylimit[7][i] = vec[i];
534*59599516SKenneth E. Jansen       }
535*59599516SKenneth E. Jansen       vec.erase(vec.begin(),vec.end());
536*59599516SKenneth E. Jansen     }
537*59599516SKenneth E. Jansen 
538*59599516SKenneth E. Jansen     if(solscalr>3){
539*59599516SKenneth E. Jansen       inpdat.epstol[5]=inp.GetValue("Scalar 4 Solver Tolerance");
540*59599516SKenneth E. Jansen       inpdat.LHSupd[5]=inp.GetValue("Number of Solves of Scalar 4 per Left-hand-side Formation");
541*59599516SKenneth E. Jansen 
542*59599516SKenneth E. Jansen       vec = inp.GetValue("Limit Scalar 4");
543*59599516SKenneth E. Jansen       for(i=0; i<3 ; i++){
544*59599516SKenneth E. Jansen         turbvar.ylimit[8][i] = vec[i];
545*59599516SKenneth E. Jansen       }
546*59599516SKenneth E. Jansen       vec.erase(vec.begin(),vec.end());
547*59599516SKenneth E. Jansen     }
548*59599516SKenneth E. Jansen 
549*59599516SKenneth E. Jansen     // DISCRETIZATION CONTROL
550*59599516SKenneth E. Jansen 
551*59599516SKenneth E. Jansen     genpar.ipord = inp.GetValue("Basis Function Order");
552*59599516SKenneth E. Jansen     if((string)inp.GetValue("Time Integration Rule") == "First Order")
553*59599516SKenneth E. Jansen       inpdat.rhoinf[0] = -1 ;
554*59599516SKenneth E. Jansen     else inpdat.rhoinf[0] = (double)inp.GetValue("Time Integration Rho Infinity");
555*59599516SKenneth E. Jansen     if((string)inp.GetValue("Predictor at Start of Step")=="Same Velocity")
556*59599516SKenneth E. Jansen       genpar.ipred = 1;
557*59599516SKenneth E. Jansen     if((string)inp.GetValue("Predictor at Start of Step")=="Zero Acceleration")
558*59599516SKenneth E. Jansen       genpar.ipred = 2;
559*59599516SKenneth E. Jansen     if((string)inp.GetValue("Predictor at Start of Step")=="Same Acceleration")
560*59599516SKenneth E. Jansen       genpar.ipred = 3;
561*59599516SKenneth E. Jansen     if((string)inp.GetValue("Predictor at Start of Step")=="Same Delta")
562*59599516SKenneth E. Jansen       genpar.ipred = 4;
563*59599516SKenneth E. Jansen 
564*59599516SKenneth E. Jansen     if((string)inp.GetValue("Weak Form") == "Galerkin")
565*59599516SKenneth E. Jansen       solpar.ivart = 1;
566*59599516SKenneth E. Jansen     if((string)inp.GetValue("Weak Form") == "SUPG")
567*59599516SKenneth E. Jansen       solpar.ivart = 2;
568*59599516SKenneth E. Jansen 
569*59599516SKenneth E. Jansen     if((string)inp.GetValue("Flow Advection Form") == "Convective")
570*59599516SKenneth E. Jansen       solpar.iconvflow = 2;
571*59599516SKenneth E. Jansen     else if((string)inp.GetValue("Flow Advection Form") == "Conservative")
572*59599516SKenneth E. Jansen       solpar.iconvflow = 1;
573*59599516SKenneth E. Jansen     if((string)inp.GetValue("Scalar Advection Form") == "Convective")
574*59599516SKenneth E. Jansen       solpar.iconvsclr = 2;
575*59599516SKenneth E. Jansen     else if((string)inp.GetValue("Scalar Advection Form") == "Conservative")
576*59599516SKenneth E. Jansen       solpar.iconvsclr = 1;
577*59599516SKenneth E. Jansen     if((string)inp.GetValue("Use Conservative Scalar Convection Velocity") == "True")
578*59599516SKenneth E. Jansen       sclrs.consrv_sclr_conv_vel = 1;
579*59599516SKenneth E. Jansen     else if((string)inp.GetValue("Use Conservative Scalar Convection Velocity") == "False")
580*59599516SKenneth E. Jansen       sclrs.consrv_sclr_conv_vel = 0;
581*59599516SKenneth E. Jansen     // TAU INPUT
582*59599516SKenneth E. Jansen     if((string)inp.GetValue("Tau Matrix") == "Diagonal-Shakib")
583*59599516SKenneth E. Jansen       genpar.itau = 0;
584*59599516SKenneth E. Jansen     else  if((string)inp.GetValue("Tau Matrix") == "Diagonal-Franca")
585*59599516SKenneth E. Jansen       genpar.itau =1;
586*59599516SKenneth E. Jansen     else if((string)inp.GetValue("Tau Matrix") == "Diagonal-Jansen(dev)")
587*59599516SKenneth E. Jansen       genpar.itau = 2;
588*59599516SKenneth E. Jansen     else if((string)inp.GetValue("Tau Matrix") == "Diagonal-Compressible")
589*59599516SKenneth E. Jansen       genpar.itau = 3;
590*59599516SKenneth E. Jansen     else if((string)inp.GetValue("Tau Matrix") == "Matrix-Mallet")
591*59599516SKenneth E. Jansen       genpar.itau = 10;
592*59599516SKenneth E. Jansen     else if((string)inp.GetValue("Tau Matrix") == "Matrix-Modal")
593*59599516SKenneth E. Jansen       genpar.itau = 11;
594*59599516SKenneth E. Jansen 
595*59599516SKenneth E. Jansen     genpar.dtsfct = inp.GetValue("Tau Time Constant");
596*59599516SKenneth E. Jansen     genpar.taucfct = inp.GetValue("Tau C Scale Factor");
597*59599516SKenneth E. Jansen 
598*59599516SKenneth E. Jansen     // FLOW DISCONTINUITY CAPTURING
599*59599516SKenneth E. Jansen 
600*59599516SKenneth E. Jansen       if((string)inp.GetValue("Discontinuity Capturing") == "Off") solpar.iDC = 0;
601*59599516SKenneth E. Jansen     else if((string)inp.GetValue("Discontinuity Capturing") == "DC-mallet") solpar.iDC = 1;
602*59599516SKenneth E. Jansen     else if((string)inp.GetValue("Discontinuity Capturing") == "DC-quadratic") solpar.iDC = 2;
603*59599516SKenneth E. Jansen    else if((string)inp.GetValue("Discontinuity Capturing") == "DC-minimum") solpar.iDC = 3;
604*59599516SKenneth E. Jansen     else {
605*59599516SKenneth E. Jansen       cout<< "Condition not defined for Discontinuity Capturing \n ";
606*59599516SKenneth E. Jansen       exit(1);
607*59599516SKenneth E. Jansen     }
608*59599516SKenneth E. Jansen 
609*59599516SKenneth E. Jansen     // SCALAR DISCONTINUITY CAPTURING
610*59599516SKenneth E. Jansen 
611*59599516SKenneth E. Jansen       vector<int> ivec = inp.GetValue("Scalar Discontinuity Capturing");
612*59599516SKenneth E. Jansen       for(i=0; i< 2; i++)  solpar.idcsclr[i] = ivec[i];
613*59599516SKenneth E. Jansen       ivec.erase(ivec.begin(),ivec.end());
614*59599516SKenneth E. Jansen 
615*59599516SKenneth E. Jansen 
616*59599516SKenneth E. Jansen //        if((string)inp.GetValue("Scalar Discontinuity Capturing") == "No") solpar.idcsclr = 0;
617*59599516SKenneth E. Jansen //      else if((string)inp.GetValue("Scalar Discontinuity Capturing") == "1") solpar.idcsclr = 1;
618*59599516SKenneth E. Jansen //   else if((string)inp.GetValue("Scalar Discontinuity Capturing") == "2") solpar.idcsclr = 2;
619*59599516SKenneth E. Jansen //   else {
620*59599516SKenneth E. Jansen //        cout<< "Condition not defined for Scalar Discontinuity Capturing \n ";
621*59599516SKenneth E. Jansen //        exit(1);
622*59599516SKenneth E. Jansen //      }
623*59599516SKenneth E. Jansen     if((string)inp.GetValue("Include Viscous Correction in Stabilization") == "True")
624*59599516SKenneth E. Jansen       {
625*59599516SKenneth E. Jansen         if(genpar.ipord == 1 ) genpar.idiff = 1;
626*59599516SKenneth E. Jansen         else genpar.idiff = 2;
627*59599516SKenneth E. Jansen       }
628*59599516SKenneth E. Jansen     else { genpar.idiff = 0;}
629*59599516SKenneth E. Jansen 
630*59599516SKenneth E. Jansen     genpar.irampViscOutlet = (int)inp.GetValue("Ramp Up Viscosity Near Outlet");
631*59599516SKenneth E. Jansen 
632*59599516SKenneth E. Jansen     genpar.istretchOutlet = (int)inp.GetValue("Stretch X Coordinate Near Outlet");
633*59599516SKenneth E. Jansen 
634*59599516SKenneth E. Jansen     genpar.iremoveStabTimeTerm = (int)inp.GetValue("Remove Time Term from Stabilization");
635*59599516SKenneth E. Jansen 
636*59599516SKenneth E. Jansen     timdat.flmpl = inp.GetValue("Lumped Mass Fraction on Left-hand-side");
637*59599516SKenneth E. Jansen     timdat.flmpr = inp.GetValue("Lumped Mass Fraction on Right-hand-side");
638*59599516SKenneth E. Jansen 
639*59599516SKenneth E. Jansen     timdat.iCFLworst = 0;
640*59599516SKenneth E. Jansen     if((string)inp.GetValue("Dump CFL") == "True")
641*59599516SKenneth E. Jansen       timdat.iCFLworst = 1;
642*59599516SKenneth E. Jansen 
643*59599516SKenneth E. Jansen     intdat.intg[0][0]=inp.GetValue("Quadrature Rule on Interior");
644*59599516SKenneth E. Jansen     intdat.intg[0][1]=inp.GetValue("Quadrature Rule on Boundary");
645*59599516SKenneth E. Jansen     genpar.ibksiz = inp.GetValue("Number of Elements Per Block");
646*59599516SKenneth E. Jansen 
647*59599516SKenneth E. Jansen     ((string)inp.GetValue("Turn Off Source Terms for Scalars")
648*59599516SKenneth E. Jansen      == "True")? sclrs.nosource=1:sclrs.nosource=0;
649*59599516SKenneth E. Jansen     sclrs.tdecay=inp.GetValue("Decay Multiplier for Scalars");
650*59599516SKenneth E. Jansen 
651*59599516SKenneth E. Jansen     // TURBULENCE MODELING PARAMETER
652*59599516SKenneth E. Jansen     int tpturb = turbvari.iles-turbvari.irans;
653*59599516SKenneth E. Jansen     int ifrule;
654*59599516SKenneth E. Jansen     if( tpturb != 0 ){
655*59599516SKenneth E. Jansen 
656*59599516SKenneth E. Jansen 
657*59599516SKenneth E. Jansen       turbvari.nohomog =inp.GetValue("Number of Homogenous Directions");
658*59599516SKenneth E. Jansen 
659*59599516SKenneth E. Jansen       if((string)inp.GetValue("Turbulence Wall Model Type") == "Slip Velocity") turbvar.itwmod = 1;
660*59599516SKenneth E. Jansen       else if((string)inp.GetValue("Turbulence Wall Model Type") == "Effective Viscosity") turbvar.itwmod = 2;
661*59599516SKenneth E. Jansen       else  turbvar.itwmod = 0;
662*59599516SKenneth E. Jansen       if (turbvari.irans < 0) turbvar.itwmod = turbvar.itwmod*(-1);
663*59599516SKenneth E. Jansen       ifrule  = inp.GetValue("Velocity Averaging Steps");
664*59599516SKenneth E. Jansen       turbvar.wtavei =(ifrule >0)? 1.0/ifrule : -1.0/ifrule;
665*59599516SKenneth E. Jansen 
666*59599516SKenneth E. Jansen       if(turbvari.iles == 1){
667*59599516SKenneth E. Jansen 
668*59599516SKenneth E. Jansen         if((string)inp.GetValue("Dynamic Model Type") == "Bardina") turbvari.iles += 10;
669*59599516SKenneth E. Jansen         else if((string)inp.GetValue("Dynamic Model Type") == "Projection") turbvari.iles += 20;
670*59599516SKenneth E. Jansen 
671*59599516SKenneth E. Jansen         ifrule = inp.GetValue("Filter Integration Rule");
672*59599516SKenneth E. Jansen         turbvari.iles += ifrule-1;
673*59599516SKenneth E. Jansen         ifrule = inp.GetValue("Dynamic Model Averaging Steps");
674*59599516SKenneth E. Jansen         turbvar.dtavei = (ifrule >0)? 1.0/ifrule : -1.0/ifrule;
675*59599516SKenneth E. Jansen         turbvar.fwr1 = inp.GetValue("Filter Width Ratio");
676*59599516SKenneth E. Jansen         turbvar.flump = inp.GetValue("Lumping Factor for Filter");
677*59599516SKenneth E. Jansen 
678*59599516SKenneth E. Jansen 
679*59599516SKenneth E. Jansen         if ((string)inp.GetValue("Model Statistics") == "True" ) {
680*59599516SKenneth E. Jansen           turbvari.modlstats = 1; }
681*59599516SKenneth E. Jansen         else {
682*59599516SKenneth E. Jansen           turbvari.modlstats = 0; }
683*59599516SKenneth E. Jansen 
684*59599516SKenneth E. Jansen         if ((string)inp.GetValue("Double Filter") == "True" ) {
685*59599516SKenneth E. Jansen           turbvari.i2filt = 1; }
686*59599516SKenneth E. Jansen         else {
687*59599516SKenneth E. Jansen           turbvari.i2filt = 0; }
688*59599516SKenneth E. Jansen 
689*59599516SKenneth E. Jansen         if ((string)inp.GetValue("Model/SUPG Dissipation") == "True" ) {
690*59599516SKenneth E. Jansen           turbvari.idis = 1; }
691*59599516SKenneth E. Jansen         else {
692*59599516SKenneth E. Jansen           turbvari.idis = 0; }
693*59599516SKenneth E. Jansen 
694*59599516SKenneth E. Jansen 
695*59599516SKenneth E. Jansen         if((string)inp.GetValue("Dynamic Model Type") == "Standard") {
696*59599516SKenneth E. Jansen 
697*59599516SKenneth E. Jansen           if((string)inp.GetValue("Dynamic Sub-Model Type") == "None")
698*59599516SKenneth E. Jansen             turbvari.isubmod = 0;
699*59599516SKenneth E. Jansen           else if((string)inp.GetValue("Dynamic Sub-Model Type") =="DFWR")
700*59599516SKenneth E. Jansen             turbvari.isubmod = 1;
701*59599516SKenneth E. Jansen           else if((string)inp.GetValue("Dynamic Sub-Model Type") =="SUPG")
702*59599516SKenneth E. Jansen             turbvari.isubmod = 2;
703*59599516SKenneth E. Jansen         }
704*59599516SKenneth E. Jansen         else if((string)inp.GetValue("Dynamic Model Type") == "Projection") {
705*59599516SKenneth E. Jansen 
706*59599516SKenneth E. Jansen           if((string)inp.GetValue("Projection Filter Type") == "Linear")
707*59599516SKenneth E. Jansen             turbvari.ifproj = 0;
708*59599516SKenneth E. Jansen           else if((string)inp.GetValue("Projection Filter Type") =="Quadratic")
709*59599516SKenneth E. Jansen             turbvari.ifproj = 1;
710*59599516SKenneth E. Jansen 
711*59599516SKenneth E. Jansen           if((string)inp.GetValue("Dynamic Sub-Model Type") == "None")
712*59599516SKenneth E. Jansen             turbvari.isubmod = 0;
713*59599516SKenneth E. Jansen           else if((string)inp.GetValue("Dynamic Sub-Model Type") =="ConsistentProj")
714*59599516SKenneth E. Jansen             turbvari.isubmod = 1;
715*59599516SKenneth E. Jansen         }
716*59599516SKenneth E. Jansen 
717*59599516SKenneth E. Jansen       }
718*59599516SKenneth E. Jansen     }
719*59599516SKenneth E. Jansen 
720*59599516SKenneth E. Jansen     // SPEBC MODELING PARAMETERS
721*59599516SKenneth E. Jansen 
722*59599516SKenneth E. Jansen     if ( (spebcvr.irscale = inp.GetValue("SPEBC Model Active")) >= 0 ){
723*59599516SKenneth E. Jansen 
724*59599516SKenneth E. Jansen       ifrule  = inp.GetValue("Velocity Averaging Steps");
725*59599516SKenneth E. Jansen       turbvar.wtavei =(ifrule >0)? 1.0/ifrule : 1.0/inpdat.nstep[0];
726*59599516SKenneth E. Jansen       spebcvr.intpres = inp.GetValue("Interpolate Pressure");
727*59599516SKenneth E. Jansen       spebcvr.plandist = inp.GetValue("Distance between Planes");
728*59599516SKenneth E. Jansen       spebcvr.thetag  = inp.GetValue("Theta Angle of Arc");
729*59599516SKenneth E. Jansen       spebcvr.ds = inp.GetValue("Distance for Velocity Averaging");
730*59599516SKenneth E. Jansen       spebcvr.tolerence = inp.GetValue("SPEBC Cylindrical Tolerance");
731*59599516SKenneth E. Jansen       spebcvr.radcyl = inp.GetValue("Radius of recycle plane");
732*59599516SKenneth E. Jansen       spebcvr.rbltin  = inp.GetValue("Inlet Boundary Layer Thickness");
733*59599516SKenneth E. Jansen       spebcvr.rvscal  = inp.GetValue("Vertical Velocity Scale Factor");
734*59599516SKenneth E. Jansen     }
735*59599516SKenneth E. Jansen 
736*59599516SKenneth E. Jansen     // CARDIOVASCULAR MODELING PARAMETERS
737*59599516SKenneth E. Jansen     if ( (string)inp.GetValue("Time Varying Boundary Conditions From File") == "True")
738*59599516SKenneth E. Jansen       nomodule.itvn = 1;
739*59599516SKenneth E. Jansen     else
740*59599516SKenneth E. Jansen       nomodule.itvn = 0;
741*59599516SKenneth E. Jansen 
742*59599516SKenneth E. Jansen     if ( nomodule.itvn ==1)
743*59599516SKenneth E. Jansen       nomodule.bcttimescale = inp.GetValue("BCT Time Scale Factor");
744*59599516SKenneth E. Jansen 
745*59599516SKenneth E. Jansen     nomodule.ipvsq=0;
746*59599516SKenneth E. Jansen     if(nomodule.icardio = inp.GetValue("Number of Coupled Surfaces")){
747*59599516SKenneth E. Jansen       if ( nomodule.icardio > MAXSURF ) {
748*59599516SKenneth E. Jansen         cout << "Number of Coupled Surfaces > MAXSURF \n";
749*59599516SKenneth E. Jansen         exit(1);
750*59599516SKenneth E. Jansen       }
751*59599516SKenneth E. Jansen       if ( (string)inp.GetValue("Pressure Coupling") == "None")
752*59599516SKenneth E. Jansen         nomodule.ipvsq=0;
753*59599516SKenneth E. Jansen       if ( (string)inp.GetValue("Pressure Coupling") == "Explicit")
754*59599516SKenneth E. Jansen         nomodule.ipvsq=1;
755*59599516SKenneth E. Jansen       if ( (string)inp.GetValue("Pressure Coupling") == "Implicit")
756*59599516SKenneth E. Jansen         nomodule.ipvsq=2;
757*59599516SKenneth E. Jansen       if ( (string)inp.GetValue("Pressure Coupling") == "P-Implicit")
758*59599516SKenneth E. Jansen         nomodule.ipvsq=3;
759*59599516SKenneth E. Jansen 
760*59599516SKenneth E. Jansen       if(nomodule.numResistSrfs=inp.GetValue("Number of Resistance Surfaces")){
761*59599516SKenneth E. Jansen           ivec = inp.GetValue("List of Resistance Surfaces");
762*59599516SKenneth E. Jansen           for(i=0;i<MAXSURF+1; i++) nomodule.nsrflistResist[i] = 0;
763*59599516SKenneth E. Jansen           for(i=0; i< nomodule.numResistSrfs; i++){
764*59599516SKenneth E. Jansen               nomodule.nsrflistResist[i+1] = ivec[i];
765*59599516SKenneth E. Jansen           }
766*59599516SKenneth E. Jansen           vec = inp.GetValue("Resistance Values");
767*59599516SKenneth E. Jansen           for(i =0; i< MAXSURF+1 ; i++) nomodule.ValueListResist[i] = 0;
768*59599516SKenneth E. Jansen           for(i =0; i< nomodule.numResistSrfs ; i++) nomodule.ValueListResist[i+1] = vec[i];
769*59599516SKenneth E. Jansen           vec.erase(vec.begin(),vec.end());
770*59599516SKenneth E. Jansen       }
771*59599516SKenneth E. Jansen       if(nomodule.numImpSrfs=inp.GetValue("Number of Impedance Surfaces")){
772*59599516SKenneth E. Jansen           ivec = inp.GetValue("List of Impedance Surfaces");
773*59599516SKenneth E. Jansen           for(i=0;i<MAXSURF+1; i++) nomodule.nsrflistImp[i] = 0;
774*59599516SKenneth E. Jansen           for(i=0; i< nomodule.numImpSrfs; i++){
775*59599516SKenneth E. Jansen               nomodule.nsrflistImp[i+1] = ivec[i];
776*59599516SKenneth E. Jansen           }
777*59599516SKenneth E. Jansen           if ( (string)inp.GetValue("Impedance From File") == "True")
778*59599516SKenneth E. Jansen               nomodule.impfile = 1; else nomodule.impfile = 0;
779*59599516SKenneth E. Jansen       }
780*59599516SKenneth E. Jansen       if(nomodule.numRCRSrfs=inp.GetValue("Number of RCR Surfaces")){
781*59599516SKenneth E. Jansen           ivec = inp.GetValue("List of RCR Surfaces");
782*59599516SKenneth E. Jansen           for(i=0;i<MAXSURF+1; i++) nomodule.nsrflistRCR[i] = 0;
783*59599516SKenneth E. Jansen           for(i=0; i< nomodule.numRCRSrfs; i++){
784*59599516SKenneth E. Jansen               nomodule.nsrflistRCR[i+1] = ivec[i];
785*59599516SKenneth E. Jansen           }
786*59599516SKenneth E. Jansen           if ( (string)inp.GetValue("RCR Values From File") == "True")
787*59599516SKenneth E. Jansen               nomodule.ircrfile = 1; else nomodule.ircrfile = 0;
788*59599516SKenneth E. Jansen       }
789*59599516SKenneth E. Jansen     }
790*59599516SKenneth E. Jansen     nomodule.ideformwall = 0;
791*59599516SKenneth E. Jansen     if((string)inp.GetValue("Deformable Wall")=="True"){
792*59599516SKenneth E. Jansen         nomodule.ideformwall = 1;
793*59599516SKenneth E. Jansen         nomodule.rhovw = inp.GetValue("Density of Vessel Wall");
794*59599516SKenneth E. Jansen         nomodule.thicknessvw = inp.GetValue("Thickness of Vessel Wall");
795*59599516SKenneth E. Jansen         nomodule.evw = inp.GetValue("Young Mod of Vessel Wall");
796*59599516SKenneth E. Jansen         nomodule.rnuvw = inp.GetValue("Poisson Ratio of Vessel Wall");
797*59599516SKenneth E. Jansen         nomodule.rshearconstantvw = inp.GetValue("Shear Constant of Vessel Wall");
798*59599516SKenneth E. Jansen         if((string)inp.GetValue("Wall Mass Matrix for LHS") == "True") nomodule.iwallmassfactor = 1;
799*59599516SKenneth E. Jansen         else nomodule.iwallmassfactor = 0;
800*59599516SKenneth E. Jansen         if((string)inp.GetValue("Wall Stiffness Matrix for LHS") == "True") nomodule.iwallstiffactor = 1;
801*59599516SKenneth E. Jansen         else nomodule.iwallstiffactor = 0;
802*59599516SKenneth E. Jansen     }
803*59599516SKenneth E. Jansen     nomodule.iviscflux = 1;
804*59599516SKenneth E. Jansen     if((string)inp.GetValue("Viscous Flux Flag") == "True") nomodule.iviscflux = 1;
805*59599516SKenneth E. Jansen     if((string)inp.GetValue("Viscous Flux Flag") == "False") nomodule.iviscflux = 0;
806*59599516SKenneth E. Jansen 
807*59599516SKenneth E. Jansen 
808*59599516SKenneth E. Jansen     // Scaling Parameters Keywords
809*59599516SKenneth E. Jansen 
810*59599516SKenneth E. Jansen     outpar.ro = inp.GetValue("Density");
811*59599516SKenneth E. Jansen     outpar.vel = inp.GetValue("Velocity");
812*59599516SKenneth E. Jansen     outpar.press = inp.GetValue("Pressure");
813*59599516SKenneth E. Jansen     outpar.temper = inp.GetValue("Temperature");
814*59599516SKenneth E. Jansen     outpar.entrop = inp.GetValue("Entropy");
815*59599516SKenneth E. Jansen 
816*59599516SKenneth E. Jansen     // Step Sequencing
817*59599516SKenneth E. Jansen 
818*59599516SKenneth E. Jansen 
819*59599516SKenneth E. Jansen     ivec = inp.GetValue("Step Construction");
820*59599516SKenneth E. Jansen     sequence.seqsize = ivec.size();
821*59599516SKenneth E. Jansen     if( sequence.seqsize > 100 || sequence.seqsize < 2 )
822*59599516SKenneth E. Jansen      cerr<<"Sequence size must be between 2 and 100 "<<endl;
823*59599516SKenneth E. Jansen 
824*59599516SKenneth E. Jansen     for(i=0; i< sequence.seqsize; i++)
825*59599516SKenneth E. Jansen       sequence.stepseq[i] = ivec[i];
826*59599516SKenneth E. Jansen 
827*59599516SKenneth E. Jansen   }
828*59599516SKenneth E. Jansen   catch ( exception &e ) {
829*59599516SKenneth E. Jansen     cout << endl << "Input exception: " << e.what() << endl << endl;
830*59599516SKenneth E. Jansen     ierr = 001;
831*59599516SKenneth E. Jansen     print_error_code(ierr);
832*59599516SKenneth E. Jansen     return ierr;
833*59599516SKenneth E. Jansen   }
834*59599516SKenneth E. Jansen 
835*59599516SKenneth E. Jansen   return ierr;
836*59599516SKenneth E. Jansen 
837*59599516SKenneth E. Jansen }
838*59599516SKenneth E. Jansen 
839*59599516SKenneth E. Jansen void print_error_code(int ierr) {
840*59599516SKenneth E. Jansen   /*
841*59599516SKenneth E. Jansen     Return Error codes:
842*59599516SKenneth E. Jansen     0xx         Input error
843*59599516SKenneth E. Jansen     1xx         Solution Control
844*59599516SKenneth E. Jansen     105         Turbulence Model not supported
845*59599516SKenneth E. Jansen 
846*59599516SKenneth E. Jansen     2xx         Material Properties
847*59599516SKenneth E. Jansen 
848*59599516SKenneth E. Jansen     3xx         Output Control
849*59599516SKenneth E. Jansen 
850*59599516SKenneth E. Jansen     4xx         Discretization Control
851*59599516SKenneth E. Jansen 
852*59599516SKenneth E. Jansen     5xx         Scaling Parameters
853*59599516SKenneth E. Jansen 
854*59599516SKenneth E. Jansen     6xx         Linear Solver Control
855*59599516SKenneth E. Jansen     601         linear solver type not supported
856*59599516SKenneth E. Jansen   */
857*59599516SKenneth E. Jansen   cout << endl << endl << "Input error detected: " << endl << endl;
858*59599516SKenneth E. Jansen   if ( ierr == 001 ) {
859*59599516SKenneth E. Jansen     cout << endl << "Input Directive not understood" << endl << endl;
860*59599516SKenneth E. Jansen   }
861*59599516SKenneth E. Jansen   if ( ierr == 105 ) {
862*59599516SKenneth E. Jansen     cout << endl << "Turbulence Model Not Supported" << endl << endl;
863*59599516SKenneth E. Jansen   }
864*59599516SKenneth E. Jansen   if ( ierr == 601 ) {
865*59599516SKenneth E. Jansen     cout << endl << "Linear Solver Type Not Supported" << endl << endl;
866*59599516SKenneth E. Jansen   }
867*59599516SKenneth E. Jansen 
868*59599516SKenneth E. Jansen }
869