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