xref: /phasta/phSolver/common/input_fform.cc (revision fcf561c152dd7c9626a95a4761c186c894f98473)
159599516SKenneth E. Jansen #include <fstream>
259599516SKenneth E. Jansen #include <stdio.h>
359599516SKenneth E. Jansen #include <stdlib.h>
459599516SKenneth E. Jansen #include <vector>
559599516SKenneth E. Jansen #include <string>
659599516SKenneth E. Jansen //MR CHANGE
759599516SKenneth E. Jansen #include <cstring>
859599516SKenneth E. Jansen //MR CHANGE END
959599516SKenneth E. Jansen 
1059599516SKenneth E. Jansen #include "Input.h"
1159599516SKenneth E. Jansen #include "common_c.h"
1259599516SKenneth E. Jansen 
1359599516SKenneth E. Jansen using namespace std; //::cout;
1459599516SKenneth E. Jansen void print_error_code(int ierr);
1559599516SKenneth E. Jansen int input_fform(char inpfname[]);
1659599516SKenneth E. Jansen int SONFATH=0;
1759599516SKenneth E. Jansen extern "C" char phasta_iotype[80];
1859599516SKenneth E. Jansen 
1959599516SKenneth E. Jansen extern "C" {
2059599516SKenneth E. Jansen 
2159599516SKenneth E. Jansen   int input_fform_(char inpfname[]) {
2259599516SKenneth E. Jansen     return input_fform(inpfname);
2359599516SKenneth E. Jansen   }
2459599516SKenneth E. Jansen 
2559599516SKenneth E. Jansen }
2659599516SKenneth E. Jansen 
2759599516SKenneth E. Jansen int input_fform(char inpfname[])
2859599516SKenneth E. Jansen {
2959599516SKenneth E. Jansen 
3059599516SKenneth E. Jansen   int ierr = 0 ;
3159599516SKenneth E. Jansen   int i,j;
3259599516SKenneth E. Jansen   char* path_to_config = 0 ;
3359599516SKenneth E. Jansen   char complete_filename[256];
3459599516SKenneth E. Jansen 
3559599516SKenneth E. Jansen   try {
3659599516SKenneth E. Jansen     // get the input file stream
3759599516SKenneth E. Jansen     path_to_config = getenv("PHASTA_CONFIG");
3859599516SKenneth E. Jansen     if(path_to_config) strcpy(complete_filename, path_to_config);
3959599516SKenneth E. Jansen     else strcpy(complete_filename,".");
4059599516SKenneth E. Jansen     strcat(complete_filename, "/input.config");
4159599516SKenneth E. Jansen     if(workfc.myrank==workfc.master) {
4259599516SKenneth E. Jansen       printf("\n Complete Filename: %s \n", complete_filename);
4359599516SKenneth E. Jansen       printf("\n Local Config: %s \n\n", inpfname);
4459599516SKenneth E. Jansen     }
4559599516SKenneth E. Jansen     string def(complete_filename);
4659599516SKenneth E. Jansen     Input inp(inpfname,def);
4759599516SKenneth E. Jansen 
4859599516SKenneth E. Jansen #ifdef AMG
4959599516SKenneth E. Jansen 
5059599516SKenneth E. Jansen 	// AMG PARAMETERS
5159599516SKenneth E. Jansen 
5259599516SKenneth E. Jansen 	if ((string)inp.GetValue("Employ AMG") == "True" ) {
5359599516SKenneth E. Jansen 
5459599516SKenneth E. Jansen 	    amgvari.irun_amg = 1;
5559599516SKenneth E. Jansen 
5659599516SKenneth E. Jansen         amgvari.irun_amg_prec = inp.GetValue("Run AMG As CG-preconditioner");
5759599516SKenneth E. Jansen 
5859599516SKenneth E. Jansen         amgvarr.strong_eps       = inp.GetValue("Strong Criterion Eps");
5959599516SKenneth E. Jansen 
6059599516SKenneth E. Jansen         amgvarr.ramg_eps         = inp.GetValue("AMG Convergence Eps");
6159599516SKenneth E. Jansen 
6259599516SKenneth E. Jansen         amgvarr.ramg_relax       = inp.GetValue("AMG Relaxation Omega");
6359599516SKenneth E. Jansen         amgvarr.ramg_trunc       = inp.GetValue("AMG Truncation Set");
6459599516SKenneth E. Jansen 
6559599516SKenneth E. Jansen         amgvari.iamg_verb        = inp.GetValue("AMG Verbosity");
6659599516SKenneth E. Jansen 
6759599516SKenneth E. Jansen         amgvari.iamg_neg_sten    = inp.GetValue("AMG Neg_Sten");
6859599516SKenneth E. Jansen 
6959599516SKenneth E. Jansen         amgvari.iamg_nlevel      = inp.GetValue("AMG Nlevel");
7059599516SKenneth E. Jansen 
7159599516SKenneth E. Jansen         amgvari.iamg_c_solver  = inp.GetValue("AMG Coarsest Solver");
7259599516SKenneth E. Jansen 
7359599516SKenneth E. Jansen         amgvari.iamg_init = 0;
7459599516SKenneth E. Jansen         amgvari.iamg_setup_frez = inp.GetValue("AMG Freeze Setup");
7559599516SKenneth E. Jansen         if ((string)inp.GetValue("AMG Interpolation Type")=="Standard")
7659599516SKenneth E. Jansen             amgvari.iamg_interp = 1;
7759599516SKenneth E. Jansen         else
7859599516SKenneth E. Jansen             amgvari.iamg_interp = 0;
7959599516SKenneth E. Jansen         amgvari.maxnev = inp.GetValue("AMG GGB nev");
8059599516SKenneth E. Jansen         amgvari.maxncv = inp.GetValue("AMG GGB ncv");
8159599516SKenneth E. Jansen 
8259599516SKenneth E. Jansen         if ((string)inp.GetValue("AMG Smoother Type")=="Gauss-Seidel")
8359599516SKenneth E. Jansen             amgvari.iamg_smoother = 1;
8459599516SKenneth E. Jansen         else if ((string)inp.GetValue("AMG Smoother Type")=="ChebyShev")
8559599516SKenneth E. Jansen             amgvari.iamg_smoother = 2;
8659599516SKenneth E. Jansen         else if ((string)inp.GetValue("AMG Smoother Type")=="MLS")
8759599516SKenneth E. Jansen             amgvari.iamg_smoother = 3;
8859599516SKenneth E. Jansen         amgvarr.ramg_chebyratio       = inp.GetValue("AMG Chebyshev Eigenvalue ratio");
8959599516SKenneth E. Jansen 
9059599516SKenneth E. Jansen         amgvari.mlsdeg = inp.GetValue("AMG MLS Degree");
9159599516SKenneth E. Jansen         amgvari.iamg_reduce = inp.GetValue("AMG Run Reduced Serial");
9259599516SKenneth E. Jansen 	}
9359599516SKenneth E. Jansen #endif
9459599516SKenneth E. Jansen 
9559599516SKenneth E. Jansen /////////////////////////////chen Sep 25 2009  Flow Control Parameters ////////
9659599516SKenneth E. Jansen // Take BC from IC at inlet
9759599516SKenneth E. Jansen       ctrlvari.iI2Binlet = (int)inp.GetValue("Take BC from IC at Inlet");
9859599516SKenneth E. Jansen       if(ctrlvari.iI2Binlet > 0 )
9959599516SKenneth E. Jansen         {
10059599516SKenneth E. Jansen          ctrlvar.inletVelX = inp.GetValue("Inlet Bulk x Velocity");
10159599516SKenneth E. Jansen         }
10259599516SKenneth E. Jansen // set uniform outlet pressure
10359599516SKenneth E. Jansen       ctrlvari.isetOutPres = (int)inp.GetValue("Set Outlet Pressure");
10459599516SKenneth E. Jansen       if(ctrlvari.isetOutPres > 0 )
10559599516SKenneth E. Jansen         {
10659599516SKenneth E. Jansen          ctrlvar.outPres1 = inp.GetValue("Uniform Outlet Pressure");
10759599516SKenneth E. Jansen         }
10859599516SKenneth E. Jansen // set initial condition
10959599516SKenneth E. Jansen       ctrlvari.isetInitial=(int)inp.GetValue("Specify Initial Conditions");
11059599516SKenneth E. Jansen       if(ctrlvari.isetInitial==1)
11159599516SKenneth E. Jansen         {
11259599516SKenneth E. Jansen          ctrlvar.xvel_ini=inp.GetValue("Initial X Velocity");
11359599516SKenneth E. Jansen          ctrlvar.yvel_ini=inp.GetValue("Initial Y Velocity");
11459599516SKenneth E. Jansen          ctrlvar.zvel_ini=inp.GetValue("Initial Z Velocity");
11559599516SKenneth E. Jansen          ctrlvar.temp_ini=inp.GetValue("Initial Temp");
11659599516SKenneth E. Jansen          ctrlvar.pres_ini=inp.GetValue("Initial Pressure");
11759599516SKenneth E. Jansen          ctrlvar.evis_ini=inp.GetValue("Initial Scalar 1");
11859599516SKenneth E. Jansen         }
11959599516SKenneth E. Jansen 
12059599516SKenneth E. Jansen ///////////////////////////////////////////////////////////////////////////////
12159599516SKenneth E. Jansen 
12259599516SKenneth E. Jansen 
12359599516SKenneth E. Jansen     // Disabled Features
12459599516SKenneth E. Jansen 
12559599516SKenneth E. Jansen     conpar.iALE = inp.GetValue("iALE");
12659599516SKenneth E. Jansen     conpar.icoord = inp.GetValue("icoord");
12759599516SKenneth E. Jansen     conpar.irs = inp.GetValue("irs");
12859599516SKenneth E. Jansen     conpar.iexec = inp.GetValue("iexec");
12959599516SKenneth E. Jansen     timpar.ntseq = inp.GetValue("ntseq");
13059599516SKenneth E. Jansen     solpar.imap = inp.GetValue("imap");
13159599516SKenneth E. Jansen 
13259599516SKenneth E. Jansen 
13359599516SKenneth E. Jansen     // Solution Control Keywords
13459599516SKenneth E. Jansen 
13559599516SKenneth E. Jansen     if((string)inp.GetValue("Equation of State") == "Incompressible") matdat.matflg[0][0] =-1 ;
13659599516SKenneth E. Jansen     if((string)inp.GetValue("Equation of State") == "Compressible") matdat.matflg[0][0] =0;
13759599516SKenneth E. Jansen     inpdat.Delt[0] = inp.GetValue("Time Step Size");
13859599516SKenneth E. Jansen     inpdat.nstep[0] = inp.GetValue("Number of Timesteps");
13959599516SKenneth E. Jansen     if((string)inp.GetValue("Viscous Control")=="Viscous") conpar.navier=1 ; else conpar.navier=0;
14059599516SKenneth E. Jansen 
14159599516SKenneth E. Jansen     if ((string)inp.GetValue("Turbulence Model") == "No-Model" ) {
14259599516SKenneth E. Jansen       turbvari.irans = 0;
14359599516SKenneth E. Jansen       turbvari.iles  = 0;
14459599516SKenneth E. Jansen     } else if ((string)inp.GetValue("Turbulence Model") == "LES" ) {
14559599516SKenneth E. Jansen       turbvari.iles  = 1;
14659599516SKenneth E. Jansen       turbvari.irans = 0;
14759599516SKenneth E. Jansen     } else if ((string)inp.GetValue("Turbulence Model") == "RANS-SA" ) {
14859599516SKenneth E. Jansen       turbvari.iles  = 0;
14959599516SKenneth E. Jansen       turbvari.irans = -1;
15059599516SKenneth E. Jansen     } else if ((string)inp.GetValue("Turbulence Model") == "RANS" ) {
15159599516SKenneth E. Jansen       turbvari.iles  = 0;
15259599516SKenneth E. Jansen       turbvari.irans = -1; // assume S-A for backward compatibility
15359599516SKenneth E. Jansen     } else if ((string)inp.GetValue("Turbulence Model") == "RANS-KE" ) {
15459599516SKenneth E. Jansen       turbvari.iles  = 0;
15559599516SKenneth E. Jansen       turbvari.irans = -2;
15659599516SKenneth E. Jansen     } else if ((string)inp.GetValue("Turbulence Model") == "DES97" ) {
15759599516SKenneth E. Jansen       turbvari.iles  = -1;
15859599516SKenneth E. Jansen       turbvari.irans = -1;
15959599516SKenneth E. Jansen     } else if ((string)inp.GetValue("Turbulence Model") == "DDES" ) {
16059599516SKenneth E. Jansen       turbvari.iles  = -2;
16159599516SKenneth E. Jansen       turbvari.irans = -1;
16259599516SKenneth E. Jansen 
16359599516SKenneth E. Jansen     } else {
16459599516SKenneth E. Jansen       cout << " Turbulence Model: Only Legal Values ( No-Model, LES, RANS-SA, RANS-KE, DES97, DDES )";
16559599516SKenneth E. Jansen       cout << endl;
16659599516SKenneth E. Jansen       exit(1);
16759599516SKenneth E. Jansen     }
16859599516SKenneth E. Jansen 
16959599516SKenneth E. Jansen  //   if (turbvari.iles*turbvari.irans!=0) turbvar.eles=
17059599516SKenneth E. Jansen  //                                          inp.GetValue("DES Edge Length");
17159599516SKenneth E. Jansen 
17259599516SKenneth E. Jansen     if (turbvari.irans<0 && turbvari.iles<0)
17359599516SKenneth E. Jansen       turbvar.DES_SA_hmin=(double)inp.GetValue("DES SA Minimum Edge Length");
17459599516SKenneth E. Jansen 
17559599516SKenneth E. Jansen     int solflow, solheat , solscalr, ilset;
17659599516SKenneth E. Jansen     ((string)inp.GetValue("Solve Flow") == "True")? solflow=1:solflow=0;
17759599516SKenneth E. Jansen     ((string)inp.GetValue("Solve Heat") == "True")? solheat=1:solheat=0;
17859599516SKenneth E. Jansen     //for compressible solheat= False so
17959599516SKenneth E. Jansen     if((string)inp.GetValue("Equation of State") == "Compressible") solheat=0;
18059599516SKenneth E. Jansen     ilset = (int)inp.GetValue("Solve Level Set");
18159599516SKenneth E. Jansen     solscalr = (int)inp.GetValue("Solve Scalars");
18259599516SKenneth E. Jansen     solscalr += ilset;
18359599516SKenneth E. Jansen     if(turbvari.irans == -1) solscalr++;
18459599516SKenneth E. Jansen     if(turbvari.irans == -2) solscalr=solscalr+2;
18559599516SKenneth E. Jansen     if ( solscalr > 4 ) {
18659599516SKenneth E. Jansen       cout << " Only Four Scalars are supported \n";
18759599516SKenneth E. Jansen       cout <<" Please reduce number of scalars \n";
18859599516SKenneth E. Jansen       exit(1);
18959599516SKenneth E. Jansen     }
19059599516SKenneth E. Jansen     inpdat.impl[0] = 10*solflow+solscalr*100+solheat;
19159599516SKenneth E. Jansen 
19259599516SKenneth E. Jansen     levlset.iLSet = ilset;
19359599516SKenneth E. Jansen     if( ilset > 0) {
19459599516SKenneth E. Jansen     levlset.epsilon_ls = inp.GetValue("Number of Elements Across Interface");
19559599516SKenneth E. Jansen     levlset.epsilon_lsd = inp.GetValue("Number of Elements Across Interface for Redistancing");
19659599516SKenneth E. Jansen     levlset.dtlset = inp.GetValue("Pseudo Time step for Redistancing");
19759599516SKenneth E. Jansen     levlset.iExpLSSclr2 = inp.GetValue("Explicit Solve for Redistance Field");
19859599516SKenneth E. Jansen     levlset.iExpLSSclr1 = inp.GetValue("Explicit Solve for Scalar 1 Field");
19959599516SKenneth E. Jansen     if ((string)inp.GetValue("Apply Volume Constraint") == "True" ) {
20059599516SKenneth E. Jansen       levlset.ivconstraint = 1; }
20159599516SKenneth E. Jansen     else if((string)inp.GetValue("Apply Volume Constraint") == "False" ) {
20259599516SKenneth E. Jansen       levlset.ivconstraint = 0; }
20359599516SKenneth E. Jansen     else {
20459599516SKenneth E. Jansen       cout << "Apply Volume Constraint: Only Legal Values (True, False) ";
20559599516SKenneth E. Jansen       cout << endl;
20659599516SKenneth E. Jansen       exit(1);
20759599516SKenneth E. Jansen     }
20859599516SKenneth E. Jansen     }
20959599516SKenneth E. Jansen 
21059599516SKenneth E. Jansen     vector<double> vec;
21159599516SKenneth E. Jansen 
21259599516SKenneth E. Jansen     // OUTPUT CONTROL KEY WORDS.
21359599516SKenneth E. Jansen 
21459599516SKenneth E. Jansen     conpar.necho = inp.GetValue("Verbosity Level");
21559599516SKenneth E. Jansen     outpar.ntout = inp.GetValue("Number of Timesteps between Restarts");
216*fcf561c1SCameron Smith     outpar.nsynciofiles = inp.GetValue("Number of SyncIO Files");
21759599516SKenneth E. Jansen     if((string)inp.GetValue("Print Statistics") == "True") outpar.ioform = 2;
21859599516SKenneth E. Jansen     else outpar.ioform = 1;
21959599516SKenneth E. Jansen 
22059599516SKenneth E. Jansen     if((string)inp.GetValue("Print Wall Fluxes") == "True") outpar.iowflux = 1;
22159599516SKenneth E. Jansen     else outpar.iowflux = 0;
22259599516SKenneth E. Jansen 
22359599516SKenneth E. Jansen     if((string)inp.GetValue("Print FieldView") == "True") outpar.iofieldv = 1;
22459599516SKenneth E. Jansen     else outpar.iofieldv = 0;
22559599516SKenneth E. Jansen 
22659599516SKenneth E. Jansen     if((string)inp.GetValue("Print ybar") == "True") outpar.ioybar = 1;
22759599516SKenneth E. Jansen     else outpar.ioybar = 0;
22859599516SKenneth E. Jansen 
22959599516SKenneth E. Jansen     if((string)inp.GetValue("Print vorticity") == "True") outpar.ivort = 1;
23059599516SKenneth E. Jansen     else outpar.ivort = 0;
23159599516SKenneth E. Jansen 
23259599516SKenneth E. Jansen     outpar.nstepsincycle = inp.GetValue("Number of Steps in a Cycle");
23359599516SKenneth E. Jansen     outpar.nphasesincycle = inp.GetValue("Number of Phases in a Cycle");
23459599516SKenneth E. Jansen     outpar.ncycles_startphaseavg = inp.GetValue("Number of Initial Cycles to Skip in Phase Average");
23559599516SKenneth E. Jansen 
23659599516SKenneth E. Jansen     strcpy( outpar.iotype , ((string)inp.GetValue("Data Block Format")).c_str());
23759599516SKenneth E. Jansen     strcpy( phasta_iotype , ((string)inp.GetValue("Data Block Format")).c_str());
23859599516SKenneth E. Jansen     SONFATH = inp.GetValue("Number of Father Nodes");
23959599516SKenneth E. Jansen 
24059599516SKenneth E. Jansen     if((string)inp.GetValue("Print Residual at End of Step") == "True") genpar.lstres = 1;
24159599516SKenneth E. Jansen     else genpar.lstres = 0;
24259599516SKenneth E. Jansen 
24359599516SKenneth E. Jansen     if((string)inp.GetValue("Print Error Indicators") == "True") turbvar.ierrcalc = 1;
24459599516SKenneth E. Jansen     else turbvar.ierrcalc = 0;
24559599516SKenneth E. Jansen 
24659599516SKenneth E. Jansen     if((string)inp.GetValue("Print Velocity Hessian") == "True") turbvar.ihessian = 1;
24759599516SKenneth E. Jansen     else turbvar.ihessian = 0;
24859599516SKenneth E. Jansen 
24959599516SKenneth E. Jansen     if ( turbvar.ierrcalc == 1 )
25059599516SKenneth E. Jansen         turbvari.ierrsmooth = inp.GetValue("Number of Error Smoothing Iterations");
25159599516SKenneth E. Jansen 
25259599516SKenneth E. Jansen     for(i=0;i<MAXSURF+1; i++) aerfrc.nsrflist[i] = 0;
25359599516SKenneth E. Jansen     int nsrfCM = inp.GetValue("Number of Force Surfaces");
25459599516SKenneth E. Jansen     if (nsrfCM > 0) {
25559599516SKenneth E. Jansen       vector<int> ivec = inp.GetValue("Surface ID's for Force Calculation");
25659599516SKenneth E. Jansen       for(i=0; i< nsrfCM; i++){
25759599516SKenneth E. Jansen         aerfrc.nsrflist[ivec[i]] = 1;
25859599516SKenneth E. Jansen         //        cout <<"surface in force list "<< ivec[i] << endl;
25959599516SKenneth E. Jansen       }
26059599516SKenneth E. Jansen       ivec.erase(ivec.begin(),ivec.end());
26159599516SKenneth E. Jansen     }
26259599516SKenneth E. Jansen 
26359599516SKenneth E. Jansen     aerfrc.isrfIM = inp.GetValue("Surface ID for Integrated Mass");
26459599516SKenneth E. Jansen 
26559599516SKenneth E. Jansen     outpar.iv_rankpercore = inp.GetValue("Ranks per core");
26659599516SKenneth E. Jansen     outpar.iv_corepernode = inp.GetValue("Cores per node");
26759599516SKenneth E. Jansen 
26859599516SKenneth E. Jansen     turbvari.iramp=0;
26959599516SKenneth E. Jansen     if((string)inp.GetValue("Ramp Inflow") == "True") turbvari.iramp=1;
27059599516SKenneth E. Jansen     if(turbvari.iramp == 1) {
27159599516SKenneth E. Jansen 	vec = inp.GetValue("Mdot Ramp Inflow Start and Stop");
27259599516SKenneth E. Jansen     	for(i=0; i<2 ; i++){
27359599516SKenneth E. Jansen         	turbvar.rampmdot[0][i]=vec[i];
27459599516SKenneth E. Jansen     	}
27559599516SKenneth E. Jansen     	vec = inp.GetValue("Mdot Ramp Lower FC Start and Stop");
27659599516SKenneth E. Jansen     	for(i=0; i<2 ; i++){
27759599516SKenneth E. Jansen          	turbvar.rampmdot[1][i]=vec[i];
27859599516SKenneth E. Jansen     	}
27959599516SKenneth E. Jansen     	vec = inp.GetValue("Mdot Ramp Upper FC Start and Stop");
28059599516SKenneth E. Jansen     	for(i=0; i<2 ; i++){
28159599516SKenneth E. Jansen         	turbvar.rampmdot[2][i]=vec[i];
28259599516SKenneth E. Jansen     	}
28359599516SKenneth E. Jansen     }
28459599516SKenneth E. Jansen 
28559599516SKenneth E. Jansen //Limiting
28659599516SKenneth E. Jansen     vec = inp.GetValue("Limit u1");
28759599516SKenneth E. Jansen     for(i=0; i<3 ; i++){
28859599516SKenneth E. Jansen       turbvar.ylimit[0][i] = vec[i];
28959599516SKenneth E. Jansen     }
29059599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
29159599516SKenneth E. Jansen 
29259599516SKenneth E. Jansen     vec = inp.GetValue("Limit u2");
29359599516SKenneth E. Jansen     for(i=0; i<3 ; i++){
29459599516SKenneth E. Jansen       turbvar.ylimit[1][i] = vec[i];
29559599516SKenneth E. Jansen     }
29659599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
29759599516SKenneth E. Jansen 
29859599516SKenneth E. Jansen     vec = inp.GetValue("Limit u3");
29959599516SKenneth E. Jansen     for(i=0; i<3 ; i++){
30059599516SKenneth E. Jansen       turbvar.ylimit[2][i] = vec[i];
30159599516SKenneth E. Jansen     }
30259599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
30359599516SKenneth E. Jansen 
30459599516SKenneth E. Jansen     vec = inp.GetValue("Limit Pressure");
30559599516SKenneth E. Jansen     for(i=0; i<3 ; i++){
30659599516SKenneth E. Jansen       turbvar.ylimit[3][i] = vec[i];
30759599516SKenneth E. Jansen     }
30859599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
30959599516SKenneth E. Jansen 
31059599516SKenneth E. Jansen     vec = inp.GetValue("Limit Temperature");
31159599516SKenneth E. Jansen     for(i=0; i<3 ; i++){
31259599516SKenneth E. Jansen       turbvar.ylimit[4][i] = vec[i];
31359599516SKenneth E. Jansen     }
31459599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
31559599516SKenneth E. Jansen 
31659599516SKenneth E. Jansen     //Material Properties Keywords
31759599516SKenneth E. Jansen     matdat.nummat = levlset.iLSet+1;
31859599516SKenneth E. Jansen     if((string)inp.GetValue("Shear Law") == "Constant Viscosity")
31959599516SKenneth E. Jansen       for(i=0; i < levlset.iLSet+1; i++) matdat.matflg[i][1] = 0;
32059599516SKenneth E. Jansen 
32159599516SKenneth E. Jansen     if((string)inp.GetValue("Bulk Viscosity Law") == "Constant Bulk Viscosity")
32259599516SKenneth E. Jansen       for(i=0; i < levlset.iLSet+1; i++) matdat.matflg[i][2] = 0;
32359599516SKenneth E. Jansen 
32459599516SKenneth E. Jansen     mmatpar.pr = inp.GetValue("Prandtl Number");
32559599516SKenneth E. Jansen 
32659599516SKenneth E. Jansen     if((string)inp.GetValue("Conductivity Law") == "Constant Conductivity")
32759599516SKenneth E. Jansen       for(i=0; i < levlset.iLSet+1; i++) matdat.matflg[i][3] = 0;
32859599516SKenneth E. Jansen 
32959599516SKenneth E. Jansen     vec = inp.GetValue("Density");
33059599516SKenneth E. Jansen     for(i=0; i< levlset.iLSet +1 ; i++){
33159599516SKenneth E. Jansen       matdat.datmat[i][0][0] = vec[i];
33259599516SKenneth E. Jansen     }
33359599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
33459599516SKenneth E. Jansen 
33559599516SKenneth E. Jansen     vec = inp.GetValue("Viscosity");
33659599516SKenneth E. Jansen     for(i=0; i< levlset.iLSet +1 ; i++){
33759599516SKenneth E. Jansen       matdat.datmat[i][1][0] = vec[i];
33859599516SKenneth E. Jansen     }
33959599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
34059599516SKenneth E. Jansen 
34159599516SKenneth E. Jansen //      vec = inp.GetValue("Specific Heat");
34259599516SKenneth E. Jansen     for(i=0; i< levlset.iLSet +1 ; i++){
34359599516SKenneth E. Jansen       matdat.datmat[i][2][0] = 0;
34459599516SKenneth E. Jansen     }
34559599516SKenneth E. Jansen //      vec.erase(vec.begin(),vec.end());
34659599516SKenneth E. Jansen 
34759599516SKenneth E. Jansen     vec = inp.GetValue("Thermal Conductivity");
34859599516SKenneth E. Jansen     for(i=0; i< levlset.iLSet +1 ; i++){
34959599516SKenneth E. Jansen       matdat.datmat[i][3][0] = vec[i];
35059599516SKenneth E. Jansen     }
35159599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
35259599516SKenneth E. Jansen 
35359599516SKenneth E. Jansen     vec = inp.GetValue("Scalar Diffusivity");
35459599516SKenneth E. Jansen     for(i=0; i< solscalr ; i++){
35559599516SKenneth E. Jansen       sclrs.scdiff[i] = vec[i];
35659599516SKenneth E. Jansen     }
35759599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
35859599516SKenneth E. Jansen 
35959599516SKenneth E. Jansen     if((string)inp.GetValue("Zero Mean Pressure") == "True")
36059599516SKenneth E. Jansen       turbvar.pzero=1;
36159599516SKenneth E. Jansen 
36259599516SKenneth E. Jansen     turbvar.rmutarget = inp.GetValue("Target Viscosity For Step NSTEP");
36359599516SKenneth E. Jansen 
36459599516SKenneth E. Jansen     if ( (string)inp.GetValue("Body Force Option") == "None" ) {
36559599516SKenneth E. Jansen       for( i=0; i< levlset.iLSet +1 ; i++)  matdat.matflg[i][4] = 0;
36659599516SKenneth E. Jansen     }
36759599516SKenneth E. Jansen     else if ( (string)inp.GetValue("Body Force Option") == "Vector" ) {
36859599516SKenneth E. Jansen       for( i=0; i< levlset.iLSet +1 ; i++)  matdat.matflg[i][4] = 1;
36959599516SKenneth E. Jansen     }
37059599516SKenneth E. Jansen     else if ( (string)inp.GetValue("Body Force Option") == "User e3source.f" ) {
37159599516SKenneth E. Jansen       for( i=0; i< levlset.iLSet +1 ; i++) matdat.matflg[i][4] = 3;
37259599516SKenneth E. Jansen     }
37359599516SKenneth E. Jansen     else if ( (string)inp.GetValue("Body Force Option") == "Boussinesq" ) {
37459599516SKenneth E. Jansen       for(i=0; i< levlset.iLSet +1 ; i++) matdat.matflg[i][4] = 2;
37559599516SKenneth E. Jansen     }
37659599516SKenneth E. Jansen     else if ( (string)inp.GetValue("Body Force Option") == "Cooling Analytic" ) {
37759599516SKenneth E. Jansen       for(i=0; i< levlset.iLSet +1 ; i++) matdat.matflg[i][4] = 4;
37859599516SKenneth E. Jansen     }
37959599516SKenneth E. Jansen     else if ( (string)inp.GetValue("Body Force Option") == "Cooling Initial Condition" ) {
38059599516SKenneth E. Jansen       for(i=0; i< levlset.iLSet +1 ; i++) matdat.matflg[i][4] = 5;
38159599516SKenneth E. Jansen     }
38259599516SKenneth E. Jansen 
38359599516SKenneth E. Jansen     // the following block of stuff is common to all cooling type sponges.
38459599516SKenneth E. Jansen     // Specific stuff belongs in the conditionals above
38559599516SKenneth E. Jansen 
38659599516SKenneth E. Jansen     if(matdat.matflg[0][4] >=4) {
38759599516SKenneth E. Jansen       spongevar.betamax = inp.GetValue("Maximum Value of Sponge Parameter");
38859599516SKenneth E. Jansen       spongevar.zinsponge = inp.GetValue("Inflow Cooling Sponge Ends at z");
38959599516SKenneth E. Jansen       spongevar.zoutsponge= inp.GetValue("Outflow Cooling Sponge Begins at z");
39059599516SKenneth E. Jansen       spongevar.radsponge = inp.GetValue("Radial Cooling Sponge Begins at r");
39159599516SKenneth E. Jansen       spongevar.grthosponge = inp.GetValue("Sponge Growth Coefficient Outflow");
39259599516SKenneth E. Jansen       spongevar.grthisponge = inp.GetValue("Sponge Growth Coefficient Inflow");
39359599516SKenneth E. Jansen 
39459599516SKenneth E. Jansen 
39559599516SKenneth E. Jansen       spongevar.spongecontinuity = 0;
39659599516SKenneth E. Jansen       spongevar.spongemomentum1 = 0;
39759599516SKenneth E. Jansen       spongevar.spongemomentum2 = 0;
39859599516SKenneth E. Jansen       spongevar.spongemomentum3 = 0;
39959599516SKenneth E. Jansen       spongevar.spongeenergy = 0;
40059599516SKenneth E. Jansen 
40159599516SKenneth E. Jansen       if((string)inp.GetValue("Sponge for Continuity Equation") == "True")
40259599516SKenneth E. Jansen 	spongevar.spongecontinuity = 1;
40359599516SKenneth E. Jansen       if((string)inp.GetValue("Sponge for x Momentum Equation") == "True")
40459599516SKenneth E. Jansen 	spongevar.spongemomentum1 = 1;
40559599516SKenneth E. Jansen       if((string)inp.GetValue("Sponge for y Momentum Equation") == "True")
40659599516SKenneth E. Jansen 	spongevar.spongemomentum2 = 1;
40759599516SKenneth E. Jansen       if((string)inp.GetValue("Sponge for z Momentum Equation") == "True")
40859599516SKenneth E. Jansen 	spongevar.spongemomentum3 = 1;
40959599516SKenneth E. Jansen       if((string)inp.GetValue("Sponge for Energy Equation") == "True")
41059599516SKenneth E. Jansen 	spongevar.spongeenergy = 1;
41159599516SKenneth E. Jansen 
41259599516SKenneth E. Jansen     }
41359599516SKenneth E. Jansen 
41459599516SKenneth E. Jansen     vec = inp.GetValue("Body Force");
41559599516SKenneth E. Jansen     for(i=0; i< levlset.iLSet +1 ; i++){
41659599516SKenneth E. Jansen       matdat.datmat[i][4][0] = vec[0+i*3];
41759599516SKenneth E. Jansen       matdat.datmat[i][4][1] = vec[1+i*3];
41859599516SKenneth E. Jansen       matdat.datmat[i][4][2] = vec[2+i*3];
41959599516SKenneth E. Jansen     }
42059599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
42159599516SKenneth E. Jansen 
42259599516SKenneth E. Jansen     vec = inp.GetValue("Body Force Pressure Gradient");
42359599516SKenneth E. Jansen     for(i=0; i< levlset.iLSet +1 ; i++){
42459599516SKenneth E. Jansen       matdat.datmat[i][6][0] = vec[0+i*3];
42559599516SKenneth E. Jansen       matdat.datmat[i][6][1] = vec[1+i*3];
42659599516SKenneth E. Jansen       matdat.datmat[i][6][2] = vec[2+i*3];
42759599516SKenneth E. Jansen     }
42859599516SKenneth E. Jansen     vec.erase(vec.begin(),vec.end());
42959599516SKenneth E. Jansen 
43059599516SKenneth E. Jansen     if ( (string)inp.GetValue("Surface Tension Option") == "No" ){
43159599516SKenneth E. Jansen         genpar.isurf = 0;
43259599516SKenneth E. Jansen     }
43359599516SKenneth E. Jansen     else if ((string)inp.GetValue("Surface Tension Option") == "Yes" ){
43459599516SKenneth E. Jansen         genpar.isurf = 1;
43559599516SKenneth E. Jansen     }
43659599516SKenneth E. Jansen     else {
43759599516SKenneth E. Jansen       cout << " Surface Tension: Only Legal Values (Yes, No) ";
43859599516SKenneth E. Jansen       cout << endl;
43959599516SKenneth E. Jansen       exit(1);
44059599516SKenneth E. Jansen     }
44159599516SKenneth E. Jansen     if( genpar.isurf > 0) {
44259599516SKenneth E. Jansen       genpar.Bo = inp.GetValue("Bond Number");
44359599516SKenneth E. Jansen     }
44459599516SKenneth E. Jansen 
44559599516SKenneth E. Jansen     genpar.EntropyPressure = inp.GetValue("Entropy Form of Pressure Constraint on Weight Space");
44659599516SKenneth E. Jansen 
44759599516SKenneth E. Jansen 
44859599516SKenneth E. Jansen     if ( (string)inp.GetValue("Rotating Frame of Reference") == "True" ) {
44959599516SKenneth E. Jansen       matdat.matflg[0][5] = 1;
45059599516SKenneth E. Jansen       vec = inp.GetValue("Rotating Frame of Reference Rotation Rate");
45159599516SKenneth E. Jansen       matdat.datmat[0][5][0] = vec[0];
45259599516SKenneth E. Jansen       matdat.datmat[0][5][1] = vec[1];
45359599516SKenneth E. Jansen       matdat.datmat[0][5][2] = vec[2];
45459599516SKenneth E. Jansen       vec.erase(vec.begin(),vec.end());
45559599516SKenneth E. Jansen     }
45659599516SKenneth E. Jansen     else {
45759599516SKenneth E. Jansen       matdat.matflg[0][5] = 0;
45859599516SKenneth E. Jansen       matdat.datmat[0][5][0] = 0.;
45959599516SKenneth E. Jansen       matdat.datmat[0][5][1] = 0.;
46059599516SKenneth E. Jansen       matdat.datmat[0][5][2] = 0.;
46159599516SKenneth E. Jansen     }
46259599516SKenneth E. Jansen 
46359599516SKenneth E. Jansen 
46459599516SKenneth E. Jansen     //Linear Solver parameters
46559599516SKenneth E. Jansen     if( (string)inp.GetValue("Solver Type") =="ACUSIM with P Projection" ){
46659599516SKenneth E. Jansen       incomp.iprjFlag = 0; incomp.ipresPrjFlag=1;}
46759599516SKenneth E. Jansen     else if ( (string)inp.GetValue("Solver Type") =="ACUSIM" ){
46859599516SKenneth E. Jansen       incomp.iprjFlag = 0; incomp.ipresPrjFlag=0;}
46959599516SKenneth E. Jansen     else if( (string)inp.GetValue("Solver Type") =="ACUSIM with Velocity Projection" ){
47059599516SKenneth E. Jansen       incomp.iprjFlag = 1; incomp.ipresPrjFlag=0;}
47159599516SKenneth E. Jansen     else if( (string)inp.GetValue("Solver Type") =="ACUSIM with Full Projection" ){
47259599516SKenneth E. Jansen       incomp.iprjFlag = 1; incomp.ipresPrjFlag=1;}
47359599516SKenneth E. Jansen     else if( (string)inp.GetValue("Solver Type") =="GMRES Matrix Free"){
47459599516SKenneth E. Jansen       inpdat.impl[0] += 10*solflow;}
47559599516SKenneth E. Jansen     else if( (string)inp.GetValue("Solver Type") =="GMRES EBE"){
47659599516SKenneth E. Jansen       inpdat.impl[0] += 20*solflow;}
47759599516SKenneth E. Jansen     //GMRES sparse is assumed default and has the value of 10, MFG 20,
47859599516SKenneth E. Jansen     // EBE 30
47959599516SKenneth E. Jansen 
48059599516SKenneth E. Jansen 
48159599516SKenneth E. Jansen     //    inpdat.niter[0] = inp.GetValue("Number of Solves per Time Step");
48259599516SKenneth E. Jansen     solpar.nGMRES = inp.GetValue("Number of GMRES Sweeps per Solve");
48359599516SKenneth E. Jansen     solpar.Kspace = inp.GetValue("Number of Krylov Vectors per GMRES Sweep");
48459599516SKenneth E. Jansen     inpdat.LHSupd[0] = inp.GetValue("Number of Solves per Left-hand-side Formation");
48559599516SKenneth E. Jansen     inpdat.epstol[0] = inp.GetValue("Tolerance on Momentum Equations");
48659599516SKenneth E. Jansen     incomp.prestol = inp.GetValue("Tolerance on ACUSIM Pressure Projection");
48759599516SKenneth E. Jansen     incomp.minIters = inp.GetValue("Minimum Number of Iterations per Nonlinear Iteration");
48859599516SKenneth E. Jansen     incomp.maxIters = inp.GetValue("Maximum Number of Iterations per Nonlinear Iteration");
48959599516SKenneth E. Jansen     inpdat.deltol[0][0]=inp.GetValue("Velocity Delta Ratio");
49059599516SKenneth E. Jansen     inpdat.deltol[1][0]=inp.GetValue("Pressure Delta Ratio");
49159599516SKenneth E. Jansen     incomp.nPrjs = inp.GetValue("Number of Velocity Projection Vectors");
49259599516SKenneth E. Jansen     incomp.nPresPrjs = inp.GetValue("Number of Pressure Projection Vectors");
49359599516SKenneth E. Jansen     incomp.iverbose = inp.GetValue("ACUSIM Verbosity Level");
49459599516SKenneth E. Jansen 
49559599516SKenneth E. Jansen     if(solheat==1){
49659599516SKenneth E. Jansen       inpdat.epstol[1]=inp.GetValue("Temperature Solver Tolerance");
49759599516SKenneth E. Jansen       inpdat.LHSupd[1]=inp.GetValue("Number of Solves of Temperature per Left-hand-side Formation");
49859599516SKenneth E. Jansen     }
49959599516SKenneth E. Jansen 
50059599516SKenneth E. Jansen     // The following is where you should put any inputs that are able to
50159599516SKenneth E. Jansen     // input differently for each scalar.  It is a little tedious in the code
50259599516SKenneth E. Jansen     // but it should make the solver.inp easier to understand. Note this will
50359599516SKenneth E. Jansen     // require some care with regression tests.
50459599516SKenneth E. Jansen 
50559599516SKenneth E. Jansen 
50659599516SKenneth E. Jansen     if(solscalr>0){
50759599516SKenneth E. Jansen       inpdat.epstol[2]=inp.GetValue("Scalar 1 Solver Tolerance");
50859599516SKenneth E. Jansen       inpdat.LHSupd[2]=inp.GetValue("Number of Solves of Scalar 1 per Left-hand-side Formation");
50959599516SKenneth E. Jansen 
51059599516SKenneth E. Jansen       vec = inp.GetValue("Limit Scalar 1");
51159599516SKenneth E. Jansen       for(i=0; i<3 ; i++){
51259599516SKenneth E. Jansen         turbvar.ylimit[5][i] = vec[i];
51359599516SKenneth E. Jansen       }
51459599516SKenneth E. Jansen       vec.erase(vec.begin(),vec.end());
51559599516SKenneth E. Jansen     }
51659599516SKenneth E. Jansen 
51759599516SKenneth E. Jansen     if(solscalr>1){
51859599516SKenneth E. Jansen       inpdat.epstol[3]=inp.GetValue("Scalar 2 Solver Tolerance");
51959599516SKenneth E. Jansen       inpdat.LHSupd[3]=inp.GetValue("Number of Solves of Scalar 2 per Left-hand-side Formation");
52059599516SKenneth E. Jansen 
52159599516SKenneth E. Jansen       vec = inp.GetValue("Limit Scalar 2");
52259599516SKenneth E. Jansen       for(i=0; i<3 ; i++){
52359599516SKenneth E. Jansen         turbvar.ylimit[6][i] = vec[i];
52459599516SKenneth E. Jansen       }
52559599516SKenneth E. Jansen       vec.erase(vec.begin(),vec.end());
52659599516SKenneth E. Jansen     }
52759599516SKenneth E. Jansen 
52859599516SKenneth E. Jansen     if(solscalr>2){
52959599516SKenneth E. Jansen       inpdat.epstol[4]=inp.GetValue("Scalar 3 Solver Tolerance");
53059599516SKenneth E. Jansen       inpdat.LHSupd[4]=inp.GetValue("Number of Solves of Scalar 3 per Left-hand-side Formation");
53159599516SKenneth E. Jansen 
53259599516SKenneth E. Jansen       vec = inp.GetValue("Limit Scalar 3");
53359599516SKenneth E. Jansen       for(i=0; i<3 ; i++){
53459599516SKenneth E. Jansen         turbvar.ylimit[7][i] = vec[i];
53559599516SKenneth E. Jansen       }
53659599516SKenneth E. Jansen       vec.erase(vec.begin(),vec.end());
53759599516SKenneth E. Jansen     }
53859599516SKenneth E. Jansen 
53959599516SKenneth E. Jansen     if(solscalr>3){
54059599516SKenneth E. Jansen       inpdat.epstol[5]=inp.GetValue("Scalar 4 Solver Tolerance");
54159599516SKenneth E. Jansen       inpdat.LHSupd[5]=inp.GetValue("Number of Solves of Scalar 4 per Left-hand-side Formation");
54259599516SKenneth E. Jansen 
54359599516SKenneth E. Jansen       vec = inp.GetValue("Limit Scalar 4");
54459599516SKenneth E. Jansen       for(i=0; i<3 ; i++){
54559599516SKenneth E. Jansen         turbvar.ylimit[8][i] = vec[i];
54659599516SKenneth E. Jansen       }
54759599516SKenneth E. Jansen       vec.erase(vec.begin(),vec.end());
54859599516SKenneth E. Jansen     }
54959599516SKenneth E. Jansen 
55059599516SKenneth E. Jansen     // DISCRETIZATION CONTROL
55159599516SKenneth E. Jansen 
55259599516SKenneth E. Jansen     genpar.ipord = inp.GetValue("Basis Function Order");
55359599516SKenneth E. Jansen     if((string)inp.GetValue("Time Integration Rule") == "First Order")
55459599516SKenneth E. Jansen       inpdat.rhoinf[0] = -1 ;
55559599516SKenneth E. Jansen     else inpdat.rhoinf[0] = (double)inp.GetValue("Time Integration Rho Infinity");
55659599516SKenneth E. Jansen     if((string)inp.GetValue("Predictor at Start of Step")=="Same Velocity")
55759599516SKenneth E. Jansen       genpar.ipred = 1;
55859599516SKenneth E. Jansen     if((string)inp.GetValue("Predictor at Start of Step")=="Zero Acceleration")
55959599516SKenneth E. Jansen       genpar.ipred = 2;
56059599516SKenneth E. Jansen     if((string)inp.GetValue("Predictor at Start of Step")=="Same Acceleration")
56159599516SKenneth E. Jansen       genpar.ipred = 3;
56259599516SKenneth E. Jansen     if((string)inp.GetValue("Predictor at Start of Step")=="Same Delta")
56359599516SKenneth E. Jansen       genpar.ipred = 4;
56459599516SKenneth E. Jansen 
56559599516SKenneth E. Jansen     if((string)inp.GetValue("Weak Form") == "Galerkin")
56659599516SKenneth E. Jansen       solpar.ivart = 1;
56759599516SKenneth E. Jansen     if((string)inp.GetValue("Weak Form") == "SUPG")
56859599516SKenneth E. Jansen       solpar.ivart = 2;
56959599516SKenneth E. Jansen 
57059599516SKenneth E. Jansen     if((string)inp.GetValue("Flow Advection Form") == "Convective")
57159599516SKenneth E. Jansen       solpar.iconvflow = 2;
57259599516SKenneth E. Jansen     else if((string)inp.GetValue("Flow Advection Form") == "Conservative")
57359599516SKenneth E. Jansen       solpar.iconvflow = 1;
57459599516SKenneth E. Jansen     if((string)inp.GetValue("Scalar Advection Form") == "Convective")
57559599516SKenneth E. Jansen       solpar.iconvsclr = 2;
57659599516SKenneth E. Jansen     else if((string)inp.GetValue("Scalar Advection Form") == "Conservative")
57759599516SKenneth E. Jansen       solpar.iconvsclr = 1;
57859599516SKenneth E. Jansen     if((string)inp.GetValue("Use Conservative Scalar Convection Velocity") == "True")
57959599516SKenneth E. Jansen       sclrs.consrv_sclr_conv_vel = 1;
58059599516SKenneth E. Jansen     else if((string)inp.GetValue("Use Conservative Scalar Convection Velocity") == "False")
58159599516SKenneth E. Jansen       sclrs.consrv_sclr_conv_vel = 0;
58259599516SKenneth E. Jansen     // TAU INPUT
58359599516SKenneth E. Jansen     if((string)inp.GetValue("Tau Matrix") == "Diagonal-Shakib")
58459599516SKenneth E. Jansen       genpar.itau = 0;
58559599516SKenneth E. Jansen     else  if((string)inp.GetValue("Tau Matrix") == "Diagonal-Franca")
58659599516SKenneth E. Jansen       genpar.itau =1;
58759599516SKenneth E. Jansen     else if((string)inp.GetValue("Tau Matrix") == "Diagonal-Jansen(dev)")
58859599516SKenneth E. Jansen       genpar.itau = 2;
58959599516SKenneth E. Jansen     else if((string)inp.GetValue("Tau Matrix") == "Diagonal-Compressible")
59059599516SKenneth E. Jansen       genpar.itau = 3;
59159599516SKenneth E. Jansen     else if((string)inp.GetValue("Tau Matrix") == "Matrix-Mallet")
59259599516SKenneth E. Jansen       genpar.itau = 10;
59359599516SKenneth E. Jansen     else if((string)inp.GetValue("Tau Matrix") == "Matrix-Modal")
59459599516SKenneth E. Jansen       genpar.itau = 11;
59559599516SKenneth E. Jansen 
59659599516SKenneth E. Jansen     genpar.dtsfct = inp.GetValue("Tau Time Constant");
59759599516SKenneth E. Jansen     genpar.taucfct = inp.GetValue("Tau C Scale Factor");
59859599516SKenneth E. Jansen 
59959599516SKenneth E. Jansen     // FLOW DISCONTINUITY CAPTURING
60059599516SKenneth E. Jansen 
60159599516SKenneth E. Jansen       if((string)inp.GetValue("Discontinuity Capturing") == "Off") solpar.iDC = 0;
60259599516SKenneth E. Jansen     else if((string)inp.GetValue("Discontinuity Capturing") == "DC-mallet") solpar.iDC = 1;
60359599516SKenneth E. Jansen     else if((string)inp.GetValue("Discontinuity Capturing") == "DC-quadratic") solpar.iDC = 2;
60459599516SKenneth E. Jansen    else if((string)inp.GetValue("Discontinuity Capturing") == "DC-minimum") solpar.iDC = 3;
60559599516SKenneth E. Jansen     else {
60659599516SKenneth E. Jansen       cout<< "Condition not defined for Discontinuity Capturing \n ";
60759599516SKenneth E. Jansen       exit(1);
60859599516SKenneth E. Jansen     }
60959599516SKenneth E. Jansen 
61059599516SKenneth E. Jansen     // SCALAR DISCONTINUITY CAPTURING
61159599516SKenneth E. Jansen 
61259599516SKenneth E. Jansen       vector<int> ivec = inp.GetValue("Scalar Discontinuity Capturing");
61359599516SKenneth E. Jansen       for(i=0; i< 2; i++)  solpar.idcsclr[i] = ivec[i];
61459599516SKenneth E. Jansen       ivec.erase(ivec.begin(),ivec.end());
61559599516SKenneth E. Jansen 
61659599516SKenneth E. Jansen 
61759599516SKenneth E. Jansen //        if((string)inp.GetValue("Scalar Discontinuity Capturing") == "No") solpar.idcsclr = 0;
61859599516SKenneth E. Jansen //      else if((string)inp.GetValue("Scalar Discontinuity Capturing") == "1") solpar.idcsclr = 1;
61959599516SKenneth E. Jansen //   else if((string)inp.GetValue("Scalar Discontinuity Capturing") == "2") solpar.idcsclr = 2;
62059599516SKenneth E. Jansen //   else {
62159599516SKenneth E. Jansen //        cout<< "Condition not defined for Scalar Discontinuity Capturing \n ";
62259599516SKenneth E. Jansen //        exit(1);
62359599516SKenneth E. Jansen //      }
62459599516SKenneth E. Jansen     if((string)inp.GetValue("Include Viscous Correction in Stabilization") == "True")
62559599516SKenneth E. Jansen       {
62659599516SKenneth E. Jansen         if(genpar.ipord == 1 ) genpar.idiff = 1;
62759599516SKenneth E. Jansen         else genpar.idiff = 2;
62859599516SKenneth E. Jansen       }
62959599516SKenneth E. Jansen     else { genpar.idiff = 0;}
63059599516SKenneth E. Jansen 
63159599516SKenneth E. Jansen     genpar.irampViscOutlet = (int)inp.GetValue("Ramp Up Viscosity Near Outlet");
63259599516SKenneth E. Jansen 
63359599516SKenneth E. Jansen     genpar.istretchOutlet = (int)inp.GetValue("Stretch X Coordinate Near Outlet");
63459599516SKenneth E. Jansen 
63559599516SKenneth E. Jansen     genpar.iremoveStabTimeTerm = (int)inp.GetValue("Remove Time Term from Stabilization");
63659599516SKenneth E. Jansen 
63759599516SKenneth E. Jansen     timdat.flmpl = inp.GetValue("Lumped Mass Fraction on Left-hand-side");
63859599516SKenneth E. Jansen     timdat.flmpr = inp.GetValue("Lumped Mass Fraction on Right-hand-side");
63959599516SKenneth E. Jansen 
64059599516SKenneth E. Jansen     timdat.iCFLworst = 0;
64159599516SKenneth E. Jansen     if((string)inp.GetValue("Dump CFL") == "True")
64259599516SKenneth E. Jansen       timdat.iCFLworst = 1;
64359599516SKenneth E. Jansen 
64459599516SKenneth E. Jansen     intdat.intg[0][0]=inp.GetValue("Quadrature Rule on Interior");
64559599516SKenneth E. Jansen     intdat.intg[0][1]=inp.GetValue("Quadrature Rule on Boundary");
64659599516SKenneth E. Jansen     genpar.ibksiz = inp.GetValue("Number of Elements Per Block");
64759599516SKenneth E. Jansen 
64859599516SKenneth E. Jansen     ((string)inp.GetValue("Turn Off Source Terms for Scalars")
64959599516SKenneth E. Jansen      == "True")? sclrs.nosource=1:sclrs.nosource=0;
65059599516SKenneth E. Jansen     sclrs.tdecay=inp.GetValue("Decay Multiplier for Scalars");
65159599516SKenneth E. Jansen 
65259599516SKenneth E. Jansen     // TURBULENCE MODELING PARAMETER
65359599516SKenneth E. Jansen     int tpturb = turbvari.iles-turbvari.irans;
65459599516SKenneth E. Jansen     int ifrule;
65559599516SKenneth E. Jansen     if( tpturb != 0 ){
65659599516SKenneth E. Jansen 
65759599516SKenneth E. Jansen 
65859599516SKenneth E. Jansen       turbvari.nohomog =inp.GetValue("Number of Homogenous Directions");
65959599516SKenneth E. Jansen 
66059599516SKenneth E. Jansen       if((string)inp.GetValue("Turbulence Wall Model Type") == "Slip Velocity") turbvar.itwmod = 1;
66159599516SKenneth E. Jansen       else if((string)inp.GetValue("Turbulence Wall Model Type") == "Effective Viscosity") turbvar.itwmod = 2;
66259599516SKenneth E. Jansen       else  turbvar.itwmod = 0;
66359599516SKenneth E. Jansen       if (turbvari.irans < 0) turbvar.itwmod = turbvar.itwmod*(-1);
66459599516SKenneth E. Jansen       ifrule  = inp.GetValue("Velocity Averaging Steps");
66559599516SKenneth E. Jansen       turbvar.wtavei =(ifrule >0)? 1.0/ifrule : -1.0/ifrule;
66659599516SKenneth E. Jansen 
66759599516SKenneth E. Jansen       if(turbvari.iles == 1){
66859599516SKenneth E. Jansen 
66959599516SKenneth E. Jansen         if((string)inp.GetValue("Dynamic Model Type") == "Bardina") turbvari.iles += 10;
67059599516SKenneth E. Jansen         else if((string)inp.GetValue("Dynamic Model Type") == "Projection") turbvari.iles += 20;
67159599516SKenneth E. Jansen 
67259599516SKenneth E. Jansen         ifrule = inp.GetValue("Filter Integration Rule");
67359599516SKenneth E. Jansen         turbvari.iles += ifrule-1;
67459599516SKenneth E. Jansen         ifrule = inp.GetValue("Dynamic Model Averaging Steps");
67559599516SKenneth E. Jansen         turbvar.dtavei = (ifrule >0)? 1.0/ifrule : -1.0/ifrule;
67659599516SKenneth E. Jansen         turbvar.fwr1 = inp.GetValue("Filter Width Ratio");
67759599516SKenneth E. Jansen         turbvar.flump = inp.GetValue("Lumping Factor for Filter");
67859599516SKenneth E. Jansen 
67959599516SKenneth E. Jansen 
68059599516SKenneth E. Jansen         if ((string)inp.GetValue("Model Statistics") == "True" ) {
68159599516SKenneth E. Jansen           turbvari.modlstats = 1; }
68259599516SKenneth E. Jansen         else {
68359599516SKenneth E. Jansen           turbvari.modlstats = 0; }
68459599516SKenneth E. Jansen 
68559599516SKenneth E. Jansen         if ((string)inp.GetValue("Double Filter") == "True" ) {
68659599516SKenneth E. Jansen           turbvari.i2filt = 1; }
68759599516SKenneth E. Jansen         else {
68859599516SKenneth E. Jansen           turbvari.i2filt = 0; }
68959599516SKenneth E. Jansen 
69059599516SKenneth E. Jansen         if ((string)inp.GetValue("Model/SUPG Dissipation") == "True" ) {
69159599516SKenneth E. Jansen           turbvari.idis = 1; }
69259599516SKenneth E. Jansen         else {
69359599516SKenneth E. Jansen           turbvari.idis = 0; }
69459599516SKenneth E. Jansen 
69559599516SKenneth E. Jansen 
69659599516SKenneth E. Jansen         if((string)inp.GetValue("Dynamic Model Type") == "Standard") {
69759599516SKenneth E. Jansen 
69859599516SKenneth E. Jansen           if((string)inp.GetValue("Dynamic Sub-Model Type") == "None")
69959599516SKenneth E. Jansen             turbvari.isubmod = 0;
70059599516SKenneth E. Jansen           else if((string)inp.GetValue("Dynamic Sub-Model Type") =="DFWR")
70159599516SKenneth E. Jansen             turbvari.isubmod = 1;
70259599516SKenneth E. Jansen           else if((string)inp.GetValue("Dynamic Sub-Model Type") =="SUPG")
70359599516SKenneth E. Jansen             turbvari.isubmod = 2;
70459599516SKenneth E. Jansen         }
70559599516SKenneth E. Jansen         else if((string)inp.GetValue("Dynamic Model Type") == "Projection") {
70659599516SKenneth E. Jansen 
70759599516SKenneth E. Jansen           if((string)inp.GetValue("Projection Filter Type") == "Linear")
70859599516SKenneth E. Jansen             turbvari.ifproj = 0;
70959599516SKenneth E. Jansen           else if((string)inp.GetValue("Projection Filter Type") =="Quadratic")
71059599516SKenneth E. Jansen             turbvari.ifproj = 1;
71159599516SKenneth E. Jansen 
71259599516SKenneth E. Jansen           if((string)inp.GetValue("Dynamic Sub-Model Type") == "None")
71359599516SKenneth E. Jansen             turbvari.isubmod = 0;
71459599516SKenneth E. Jansen           else if((string)inp.GetValue("Dynamic Sub-Model Type") =="ConsistentProj")
71559599516SKenneth E. Jansen             turbvari.isubmod = 1;
71659599516SKenneth E. Jansen         }
71759599516SKenneth E. Jansen 
71859599516SKenneth E. Jansen       }
71959599516SKenneth E. Jansen     }
72059599516SKenneth E. Jansen 
72159599516SKenneth E. Jansen     // SPEBC MODELING PARAMETERS
72259599516SKenneth E. Jansen 
72359599516SKenneth E. Jansen     if ( (spebcvr.irscale = inp.GetValue("SPEBC Model Active")) >= 0 ){
72459599516SKenneth E. Jansen 
72559599516SKenneth E. Jansen       ifrule  = inp.GetValue("Velocity Averaging Steps");
72659599516SKenneth E. Jansen       turbvar.wtavei =(ifrule >0)? 1.0/ifrule : 1.0/inpdat.nstep[0];
72759599516SKenneth E. Jansen       spebcvr.intpres = inp.GetValue("Interpolate Pressure");
72859599516SKenneth E. Jansen       spebcvr.plandist = inp.GetValue("Distance between Planes");
72959599516SKenneth E. Jansen       spebcvr.thetag  = inp.GetValue("Theta Angle of Arc");
73059599516SKenneth E. Jansen       spebcvr.ds = inp.GetValue("Distance for Velocity Averaging");
73159599516SKenneth E. Jansen       spebcvr.tolerence = inp.GetValue("SPEBC Cylindrical Tolerance");
73259599516SKenneth E. Jansen       spebcvr.radcyl = inp.GetValue("Radius of recycle plane");
73359599516SKenneth E. Jansen       spebcvr.rbltin  = inp.GetValue("Inlet Boundary Layer Thickness");
73459599516SKenneth E. Jansen       spebcvr.rvscal  = inp.GetValue("Vertical Velocity Scale Factor");
73559599516SKenneth E. Jansen     }
73659599516SKenneth E. Jansen 
73759599516SKenneth E. Jansen     // CARDIOVASCULAR MODELING PARAMETERS
73859599516SKenneth E. Jansen     if ( (string)inp.GetValue("Time Varying Boundary Conditions From File") == "True")
73959599516SKenneth E. Jansen       nomodule.itvn = 1;
74059599516SKenneth E. Jansen     else
74159599516SKenneth E. Jansen       nomodule.itvn = 0;
74259599516SKenneth E. Jansen 
74359599516SKenneth E. Jansen     if ( nomodule.itvn ==1)
74459599516SKenneth E. Jansen       nomodule.bcttimescale = inp.GetValue("BCT Time Scale Factor");
74559599516SKenneth E. Jansen 
74659599516SKenneth E. Jansen     nomodule.ipvsq=0;
74759599516SKenneth E. Jansen     if(nomodule.icardio = inp.GetValue("Number of Coupled Surfaces")){
74859599516SKenneth E. Jansen       if ( nomodule.icardio > MAXSURF ) {
74959599516SKenneth E. Jansen         cout << "Number of Coupled Surfaces > MAXSURF \n";
75059599516SKenneth E. Jansen         exit(1);
75159599516SKenneth E. Jansen       }
75259599516SKenneth E. Jansen       if ( (string)inp.GetValue("Pressure Coupling") == "None")
75359599516SKenneth E. Jansen         nomodule.ipvsq=0;
75459599516SKenneth E. Jansen       if ( (string)inp.GetValue("Pressure Coupling") == "Explicit")
75559599516SKenneth E. Jansen         nomodule.ipvsq=1;
75659599516SKenneth E. Jansen       if ( (string)inp.GetValue("Pressure Coupling") == "Implicit")
75759599516SKenneth E. Jansen         nomodule.ipvsq=2;
75859599516SKenneth E. Jansen       if ( (string)inp.GetValue("Pressure Coupling") == "P-Implicit")
75959599516SKenneth E. Jansen         nomodule.ipvsq=3;
76059599516SKenneth E. Jansen 
76159599516SKenneth E. Jansen       if(nomodule.numResistSrfs=inp.GetValue("Number of Resistance Surfaces")){
76259599516SKenneth E. Jansen           ivec = inp.GetValue("List of Resistance Surfaces");
76359599516SKenneth E. Jansen           for(i=0;i<MAXSURF+1; i++) nomodule.nsrflistResist[i] = 0;
76459599516SKenneth E. Jansen           for(i=0; i< nomodule.numResistSrfs; i++){
76559599516SKenneth E. Jansen               nomodule.nsrflistResist[i+1] = ivec[i];
76659599516SKenneth E. Jansen           }
76759599516SKenneth E. Jansen           vec = inp.GetValue("Resistance Values");
76859599516SKenneth E. Jansen           for(i =0; i< MAXSURF+1 ; i++) nomodule.ValueListResist[i] = 0;
76959599516SKenneth E. Jansen           for(i =0; i< nomodule.numResistSrfs ; i++) nomodule.ValueListResist[i+1] = vec[i];
77059599516SKenneth E. Jansen           vec.erase(vec.begin(),vec.end());
77159599516SKenneth E. Jansen       }
77259599516SKenneth E. Jansen       if(nomodule.numImpSrfs=inp.GetValue("Number of Impedance Surfaces")){
77359599516SKenneth E. Jansen           ivec = inp.GetValue("List of Impedance Surfaces");
77459599516SKenneth E. Jansen           for(i=0;i<MAXSURF+1; i++) nomodule.nsrflistImp[i] = 0;
77559599516SKenneth E. Jansen           for(i=0; i< nomodule.numImpSrfs; i++){
77659599516SKenneth E. Jansen               nomodule.nsrflistImp[i+1] = ivec[i];
77759599516SKenneth E. Jansen           }
77859599516SKenneth E. Jansen           if ( (string)inp.GetValue("Impedance From File") == "True")
77959599516SKenneth E. Jansen               nomodule.impfile = 1; else nomodule.impfile = 0;
78059599516SKenneth E. Jansen       }
78159599516SKenneth E. Jansen       if(nomodule.numRCRSrfs=inp.GetValue("Number of RCR Surfaces")){
78259599516SKenneth E. Jansen           ivec = inp.GetValue("List of RCR Surfaces");
78359599516SKenneth E. Jansen           for(i=0;i<MAXSURF+1; i++) nomodule.nsrflistRCR[i] = 0;
78459599516SKenneth E. Jansen           for(i=0; i< nomodule.numRCRSrfs; i++){
78559599516SKenneth E. Jansen               nomodule.nsrflistRCR[i+1] = ivec[i];
78659599516SKenneth E. Jansen           }
78759599516SKenneth E. Jansen           if ( (string)inp.GetValue("RCR Values From File") == "True")
78859599516SKenneth E. Jansen               nomodule.ircrfile = 1; else nomodule.ircrfile = 0;
78959599516SKenneth E. Jansen       }
79059599516SKenneth E. Jansen     }
79159599516SKenneth E. Jansen     nomodule.ideformwall = 0;
79259599516SKenneth E. Jansen     if((string)inp.GetValue("Deformable Wall")=="True"){
79359599516SKenneth E. Jansen         nomodule.ideformwall = 1;
79459599516SKenneth E. Jansen         nomodule.rhovw = inp.GetValue("Density of Vessel Wall");
79559599516SKenneth E. Jansen         nomodule.thicknessvw = inp.GetValue("Thickness of Vessel Wall");
79659599516SKenneth E. Jansen         nomodule.evw = inp.GetValue("Young Mod of Vessel Wall");
79759599516SKenneth E. Jansen         nomodule.rnuvw = inp.GetValue("Poisson Ratio of Vessel Wall");
79859599516SKenneth E. Jansen         nomodule.rshearconstantvw = inp.GetValue("Shear Constant of Vessel Wall");
79959599516SKenneth E. Jansen         if((string)inp.GetValue("Wall Mass Matrix for LHS") == "True") nomodule.iwallmassfactor = 1;
80059599516SKenneth E. Jansen         else nomodule.iwallmassfactor = 0;
80159599516SKenneth E. Jansen         if((string)inp.GetValue("Wall Stiffness Matrix for LHS") == "True") nomodule.iwallstiffactor = 1;
80259599516SKenneth E. Jansen         else nomodule.iwallstiffactor = 0;
80359599516SKenneth E. Jansen     }
80459599516SKenneth E. Jansen     nomodule.iviscflux = 1;
80559599516SKenneth E. Jansen     if((string)inp.GetValue("Viscous Flux Flag") == "True") nomodule.iviscflux = 1;
80659599516SKenneth E. Jansen     if((string)inp.GetValue("Viscous Flux Flag") == "False") nomodule.iviscflux = 0;
80759599516SKenneth E. Jansen 
80859599516SKenneth E. Jansen 
80959599516SKenneth E. Jansen     // Scaling Parameters Keywords
81059599516SKenneth E. Jansen 
81159599516SKenneth E. Jansen     outpar.ro = inp.GetValue("Density");
81259599516SKenneth E. Jansen     outpar.vel = inp.GetValue("Velocity");
81359599516SKenneth E. Jansen     outpar.press = inp.GetValue("Pressure");
81459599516SKenneth E. Jansen     outpar.temper = inp.GetValue("Temperature");
81559599516SKenneth E. Jansen     outpar.entrop = inp.GetValue("Entropy");
81659599516SKenneth E. Jansen 
81759599516SKenneth E. Jansen     // Step Sequencing
81859599516SKenneth E. Jansen 
81959599516SKenneth E. Jansen 
82059599516SKenneth E. Jansen     ivec = inp.GetValue("Step Construction");
82159599516SKenneth E. Jansen     sequence.seqsize = ivec.size();
82259599516SKenneth E. Jansen     if( sequence.seqsize > 100 || sequence.seqsize < 2 )
82359599516SKenneth E. Jansen      cerr<<"Sequence size must be between 2 and 100 "<<endl;
82459599516SKenneth E. Jansen 
82559599516SKenneth E. Jansen     for(i=0; i< sequence.seqsize; i++)
82659599516SKenneth E. Jansen       sequence.stepseq[i] = ivec[i];
82759599516SKenneth E. Jansen 
82859599516SKenneth E. Jansen   }
82959599516SKenneth E. Jansen   catch ( exception &e ) {
83059599516SKenneth E. Jansen     cout << endl << "Input exception: " << e.what() << endl << endl;
83159599516SKenneth E. Jansen     ierr = 001;
83259599516SKenneth E. Jansen     print_error_code(ierr);
83359599516SKenneth E. Jansen     return ierr;
83459599516SKenneth E. Jansen   }
83559599516SKenneth E. Jansen 
83659599516SKenneth E. Jansen   return ierr;
83759599516SKenneth E. Jansen 
83859599516SKenneth E. Jansen }
83959599516SKenneth E. Jansen 
84059599516SKenneth E. Jansen void print_error_code(int ierr) {
84159599516SKenneth E. Jansen   /*
84259599516SKenneth E. Jansen     Return Error codes:
84359599516SKenneth E. Jansen     0xx         Input error
84459599516SKenneth E. Jansen     1xx         Solution Control
84559599516SKenneth E. Jansen     105         Turbulence Model not supported
84659599516SKenneth E. Jansen 
84759599516SKenneth E. Jansen     2xx         Material Properties
84859599516SKenneth E. Jansen 
84959599516SKenneth E. Jansen     3xx         Output Control
85059599516SKenneth E. Jansen 
85159599516SKenneth E. Jansen     4xx         Discretization Control
85259599516SKenneth E. Jansen 
85359599516SKenneth E. Jansen     5xx         Scaling Parameters
85459599516SKenneth E. Jansen 
85559599516SKenneth E. Jansen     6xx         Linear Solver Control
85659599516SKenneth E. Jansen     601         linear solver type not supported
85759599516SKenneth E. Jansen   */
85859599516SKenneth E. Jansen   cout << endl << endl << "Input error detected: " << endl << endl;
85959599516SKenneth E. Jansen   if ( ierr == 001 ) {
86059599516SKenneth E. Jansen     cout << endl << "Input Directive not understood" << endl << endl;
86159599516SKenneth E. Jansen   }
86259599516SKenneth E. Jansen   if ( ierr == 105 ) {
86359599516SKenneth E. Jansen     cout << endl << "Turbulence Model Not Supported" << endl << endl;
86459599516SKenneth E. Jansen   }
86559599516SKenneth E. Jansen   if ( ierr == 601 ) {
86659599516SKenneth E. Jansen     cout << endl << "Linear Solver Type Not Supported" << endl << endl;
86759599516SKenneth E. Jansen   }
86859599516SKenneth E. Jansen 
86959599516SKenneth E. Jansen }
870