00001 /*---------------------------------------------------------------------------*\ 00002 ========= | 00003 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox 00004 \\ / O peration | 00005 \\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd. 00006 \\/ M anipulation | 00007 ------------------------------------------------------------------------------- 00008 License 00009 This file is part of OpenFOAM. 00010 00011 OpenFOAM is free software: you can redistribute it and/or modify it 00012 under the terms of the GNU General Public License as published by 00013 the Free Software Foundation, either version 3 of the License, or 00014 (at your option) any later version. 00015 00016 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT 00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or 00018 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License 00019 for more details. 00020 00021 You should have received a copy of the GNU General Public License 00022 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>. 00023 00024 Application 00025 shallowWaterFoam 00026 00027 Description 00028 Transient solver for inviscid shallow-water equations with rotation. 00029 00030 If the geometry is 3D then it is assumed to be one layers of cells and the 00031 component of the velocity normal to gravity is removed. 00032 00033 Usage 00034 - shallowWaterFoam [OPTION] 00035 00036 @param -case <dir> \n 00037 Specify the case directory 00038 00039 @param -parallel \n 00040 Run the case in parallel 00041 00042 @param -help \n 00043 Display short usage message 00044 00045 @param -doc \n 00046 Display Doxygen documentation page 00047 00048 @param -srcDoc \n 00049 Display source code 00050 00051 \*---------------------------------------------------------------------------*/ 00052 00053 #include <finiteVolume/fvCFD.H> 00054 00055 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00056 00057 int main(int argc, char *argv[]) 00058 { 00059 #include <OpenFOAM/setRootCase.H> 00060 #include <OpenFOAM/createTime.H> 00061 #include <OpenFOAM/createMesh.H> 00062 #include "readGravitationalAcceleration.H" 00063 #include "createFields.H" 00064 00065 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00066 00067 Info<< "\nStarting time loop\n" << endl; 00068 00069 while (runTime.loop()) 00070 { 00071 Info<< "\n Time = " << runTime.timeName() << nl << endl; 00072 00073 #include <finiteVolume/readPISOControls.H> 00074 #include "CourantNo.H" 00075 00076 for (int ucorr=0; ucorr<nOuterCorr; ucorr++) 00077 { 00078 surfaceScalarField phiv("phiv", phi/fvc::interpolate(h)); 00079 00080 fvVectorMatrix hUEqn 00081 ( 00082 fvm::ddt(hU) 00083 + fvm::div(phiv, hU) 00084 ); 00085 00086 hUEqn.relax(); 00087 00088 if (momentumPredictor) 00089 { 00090 if (rotating) 00091 { 00092 solve(hUEqn + (F ^ hU) == -magg*h*fvc::grad(h + h0)); 00093 } 00094 else 00095 { 00096 solve(hUEqn == -magg*h*fvc::grad(h + h0)); 00097 } 00098 00099 // Constrain the momentum to be in the geometry if 3D geometry 00100 if (mesh.nGeometricD() == 3) 00101 { 00102 hU -= (gHat & hU)*gHat; 00103 hU.correctBoundaryConditions(); 00104 } 00105 } 00106 00107 // --- PISO loop 00108 for (int corr=0; corr<nCorr; corr++) 00109 { 00110 surfaceScalarField hf = fvc::interpolate(h); 00111 volScalarField rUA = 1.0/hUEqn.A(); 00112 surfaceScalarField ghrUAf = magg*fvc::interpolate(h*rUA); 00113 00114 surfaceScalarField phih0 = ghrUAf*mesh.magSf()*fvc::snGrad(h0); 00115 00116 if (rotating) 00117 { 00118 hU = rUA*(hUEqn.H() - (F ^ hU)); 00119 } 00120 else 00121 { 00122 hU = rUA*hUEqn.H(); 00123 } 00124 00125 phi = (fvc::interpolate(hU) & mesh.Sf()) 00126 + fvc::ddtPhiCorr(rUA, h, hU, phi) 00127 - phih0; 00128 00129 for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++) 00130 { 00131 fvScalarMatrix hEqn 00132 ( 00133 fvm::ddt(h) 00134 + fvc::div(phi) 00135 - fvm::laplacian(ghrUAf, h) 00136 ); 00137 00138 if (ucorr < nOuterCorr-1 || corr < nCorr-1) 00139 { 00140 hEqn.solve(); 00141 } 00142 else 00143 { 00144 hEqn.solve(mesh.solver(h.name() + "Final")); 00145 } 00146 00147 if (nonOrth == nNonOrthCorr) 00148 { 00149 phi += hEqn.flux(); 00150 } 00151 } 00152 00153 hU -= rUA*h*magg*fvc::grad(h + h0); 00154 00155 // Constrain the momentum to be in the geometry if 3D geometry 00156 if (mesh.nGeometricD() == 3) 00157 { 00158 hU -= (gHat & hU)*gHat; 00159 } 00160 00161 hU.correctBoundaryConditions(); 00162 } 00163 } 00164 00165 U == hU/h; 00166 hTotal == h + h0; 00167 00168 runTime.write(); 00169 00170 Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s" 00171 << " ClockTime = " << runTime.elapsedClockTime() << " s" 00172 << nl << endl; 00173 } 00174 00175 Info<< "End\n" << endl; 00176 00177 return 0; 00178 } 00179 00180 00181 // ************************ vim: set sw=4 sts=4 et: ************************ //