1 // @HEADER
2 // ************************************************************************
3 //
4 //               Rapid Optimization Library (ROL) Package
5 //                 Copyright (2014) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact lead developers:
38 //              Drew Kouri   (dpkouri@sandia.gov) and
39 //              Denis Ridzal (dridzal@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
44 /*! \file  example_02.cpp
45     \brief Shows how to solve the binary advection-diffusion control problem
46            using sum-up rounding.
47 */
48 
49 #include "Teuchos_GlobalMPISession.hpp"
50 #include "Tpetra_Core.hpp"
51 #include "Tpetra_Version.hpp"
52 
53 #include <iostream>
54 //#include <fenv.h>
55 
56 #include "ROL_Stream.hpp"
57 #include "ROL_ParameterList.hpp"
58 #include "ROL_OptimizationSolver.hpp"
59 #include "opfactory.hpp"
60 #include "hilbert.hpp"
61 
main(int argc,char * argv[])62 int main(int argc, char *argv[]) {
63 //  feenableexcept(FE_DIVBYZERO | FE_INVALID | FE_OVERFLOW);
64   using RealT = double;
65 
66   /*** Initialize communicator. ***/
67   Teuchos::GlobalMPISession mpiSession(&argc, &argv);
68   ROL::Ptr<const Teuchos::Comm<int> > comm
69     = Tpetra::getDefaultComm();
70 
71   // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
72   const int myRank = comm->getRank();
73   ROL::Ptr<std::ostream> outStream = ROL::makeStreamPtr( std::cout, (argc > 1) && (myRank==0) );
74 
75   int errorFlag  = 0;
76 
77   // *** Example body.
78   try {
79 
80     /*** Read in XML input ***/
81     std::string filename = "input_ex01.xml";
82     ROL::Ptr<ROL::ParameterList> parlist = ROL::getParametersFromXmlFile(filename);
83     Teuchos::RCP<Teuchos::ParameterList> tparlist = Teuchos::rcp( new Teuchos::ParameterList() );
84     Teuchos::updateParametersFromXmlFile( filename, tparlist.ptr() );
85 
86     /*************************************************************************/
87     /***************** BUILD OPTIMIZATION PROBLEM ****************************/
88     /*************************************************************************/
89     ROL::Ptr<BinaryAdvDiffFactory<RealT>> factory;
90     ROL::Ptr<ROL::OptimizationProblem<RealT>> problem;
91     ROL::Ptr<ROL::OptimizationSolver<RealT>> solver;
92     ROL::Ptr<ROL::Vector<RealT>> z, omega, uomega, uz, du;
93     RealT err(0);
94     int order = parlist->sublist("Problem").get("Hilbert Curve Order",6);
95     for (int i = 0; i < order; ++i) {
96       parlist->sublist("Problem").set("Hilbert Curve Order",i+1);
97       factory = ROL::makePtr<BinaryAdvDiffFactory<RealT>>(*parlist,comm,outStream);
98       bool checkDeriv = parlist->sublist("Problem").get("Check Derivatives",false);
99       if (checkDeriv) factory->check(*outStream);
100 
101       /*************************************************************************/
102       /***************** SOLVE OPTIMIZATION PROBLEM ****************************/
103       /*************************************************************************/
104       problem = factory->build();
105       if (checkDeriv) problem->check(*outStream);
106       if (i==0) {
107         factory->getAssembler()->printMeshData(*outStream);
108       }
109       solver = ROL::makePtr<ROL::OptimizationSolver<RealT>>(*problem, *parlist);
110       solver->solve(*outStream);
111       z = problem->getSolutionVector();
112       factory->getState(uz,z);
113 
114       // Sum Up Rounding
115       omega = z->clone(); omega->zero();
116       std::vector<RealT> &zdata = *ROL::dynamicPtrCast<PDE_OptVector<RealT>>(z)->getParameter()->getVector();
117       std::vector<RealT> &odata = *ROL::dynamicPtrCast<PDE_OptVector<RealT>>(omega)->getParameter()->getVector();
118       int n = std::pow(2,i+1), d(0);
119       RealT g1(0), g2(0);
120       for (int j = 0; j < n; ++j) {
121         for (int k = 0; k < n; ++k) {
122           hilbert::xy2d(i+1,j,k,d);
123           g1 += zdata[d];
124           g2 += static_cast<RealT>(1) - zdata[d];
125           if (g1 >= g2) {
126             odata[d] = static_cast<RealT>(1);
127             g1 -= static_cast<RealT>(1);
128           }
129           else {
130             g2 -= static_cast<RealT>(1);
131           }
132         }
133       }
134 
135       // Solve PDE
136       factory->getState(uomega,omega);
137 
138       // Print
139       std::stringstream uname, zname, oname, xname, yname;
140       uname << "state_" << i+1 << ".txt";
141       zname << "control_relaxed_" << i+1 << ".txt";
142       oname << "control_binary_" << i+1 << ".txt";
143       xname << "X_" << i+1 << ".txt";
144       yname << "Y_" << i+1 << ".txt";
145       factory->getAssembler()->outputTpetraVector(ROL::dynamicPtrCast<ROL::TpetraMultiVector<RealT>>(uomega)->getVector(),uname.str());
146       std::ofstream zfile, ofile, xfile, yfile;
147       zfile.open(zname.str());
148       ofile.open(oname.str());
149       xfile.open(xname.str());
150       yfile.open(yname.str());
151       int x(0), y(0);
152       for (unsigned j = 0; j < n*n; ++j) {
153         zfile << zdata[j] << std::endl;
154         ofile << odata[j] << std::endl;
155         hilbert::d2xy(i+1, j, x, y);
156         xfile << x << std::endl;
157         yfile << y << std::endl;
158       }
159       zfile.close();
160       ofile.close();
161       xfile.close();
162       yfile.close();
163 
164       du = uz->clone();
165       du->set(*uz); du->axpy(static_cast<RealT>(-1),*uomega);
166       err = du->norm();
167       *outStream << "State Error: " << err << std::endl;
168     }
169     factory->print(*outStream);
170 
171     //*outStream << "OPTIMAL CONTROLS" << std::endl;
172     //for (int i=0; i<controlDim; ++i) {
173     //  *outStream << (*z_ptr)[i] << "  ";
174     //}
175     //*outStream << std::endl;
176   }
177   catch (std::logic_error& err) {
178     *outStream << err.what() << "\n";
179     errorFlag = -1000;
180   }; // end try
181 
182   if (errorFlag != 0)
183     std::cout << "End Result: TEST FAILED\n";
184   else
185     std::cout << "End Result: TEST PASSED\n";
186 
187   return 0;
188 }
189