/*---------------------------------------------------------------------------*\ ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗ ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║ ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║ ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝ ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝ ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝ * In real Time Highly Advanced Computational Applications for Finite Volumes * Copyright (C) 2017 by the ITHACA-FV authors ------------------------------------------------------------------------------- License This file is part of ITHACA-FV ITHACA-FV is free software: you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. ITHACA-FV is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details. You should have received a copy of the GNU Lesser General Public License along with ITHACA-FV. If not, see . Description Example of steady NS Reduction Problem solved by the use of the SIMPLE algorithm SourceFiles 12simpleSteadyNS.C \*---------------------------------------------------------------------------*/ #include "SteadyNSSimple.H" #include "ITHACAstream.H" #include "ITHACAPOD.H" #include "ReducedSimpleSteadyNS.H" #include "forces.H" #include "IOmanip.H" class tutorial12 : public SteadyNSSimple { public: /// Constructor explicit tutorial12(int argc, char* argv[]) : SteadyNSSimple(argc, argv), U(_U()), p(_p()), phi(_phi()) {} /// Velocity field volVectorField& U; /// Pressure field volScalarField& p; /// surfaceScalarField& phi; /// Perform an Offline solve void offlineSolve() { Vector inl(0, 0, 0); List mu_now(1); // if the offline solution is already performed read the fields if (offline) { ITHACAstream::read_fields(Ufield, U, "./ITHACAoutput/Offline/"); ITHACAstream::read_fields(Pfield, p, "./ITHACAoutput/Offline/"); mu_samples = ITHACAstream::readMatrix("./ITHACAoutput/Offline/mu_samples_mat.txt"); } else { Vector Uinl(1, 0, 0); label BCind = 0; for (label i = 0; i < mu.cols(); i++) { mu_now[0] = mu(0, i); change_viscosity(mu(0, i)); assignIF(U, Uinl); truthSolve2(mu_now); } } } }; int main(int argc, char* argv[]) { // Construct the tutorial object tutorial12 example(argc, argv); // Read some parameters from file ITHACAparameters* para = ITHACAparameters::getInstance(example._mesh(), example._runTime()); int NmodesUout = para->ITHACAdict->lookupOrDefault("NmodesUout", 15); int NmodesPout = para->ITHACAdict->lookupOrDefault("NmodesPout", 15); int NmodesUproj = para->ITHACAdict->lookupOrDefault("NmodesUproj", 10); int NmodesPproj = para->ITHACAdict->lookupOrDefault("NmodesPproj", 10); // Read the par file where the parameters are stored word filename("./par"); example.mu = ITHACAstream::readMatrix(filename); // Set the inlet boundaries patch 0 directions x and y example.inletIndex.resize(1, 2); example.inletIndex(0, 0) = 0; example.inletIndex(0, 1) = 0; // Perform the offline solve example.offlineSolve(); ITHACAstream::read_fields(example.liftfield, example.U, "./lift/"); // Homogenize the snapshots example.computeLift(example.Ufield, example.liftfield, example.Uomfield); // Perform POD on velocity and pressure and store the first 10 modes ITHACAPOD::getModes(example.Uomfield, example.Umodes, example._U().name(), example.podex, 0, 0, NmodesUout); ITHACAPOD::getModes(example.Pfield, example.Pmodes, example._p().name(), example.podex, 0, 0, NmodesPout); // Create the reduced object reducedSimpleSteadyNS reduced(example); PtrList U_rec_list; PtrList P_rec_list; // Reads inlet volocities boundary conditions. word vel_file(para->ITHACAdict->lookup("online_velocities")); Eigen::MatrixXd vel = ITHACAstream::readMatrix(vel_file); //Perform the online solutions for (label k = 0; k < (example.mu).size(); k++) { scalar mu_now = example.mu(0, k); example.change_viscosity(mu_now); reduced.setOnlineVelocity(vel); reduced.solveOnline_Simple(mu_now, NmodesUproj, NmodesPproj); } return 0; }