Stokhos Package Browser (Single Doxygen Collection) Version of the Day
Loading...
Searching...
No Matches
Stokhos_SGInverseModelEvaluator.cpp
Go to the documentation of this file.
1// @HEADER
2// ***********************************************************************
3//
4// Stokhos Package
5// Copyright (2009) 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 Eric T. Phipps (etphipp@sandia.gov).
38//
39// ***********************************************************************
40// @HEADER
41
47#include "Epetra_Map.h"
48#include "Teuchos_Assert.hpp"
49
51 const Teuchos::RCP<EpetraExt::ModelEvaluator>& me_,
52 const Teuchos::Array<int>& sg_p_index_map_,
53 const Teuchos::Array<int>& sg_g_index_map_,
54 const Teuchos::Array< Teuchos::RCP<const Epetra_Map> >& base_g_maps_) :
55 me(me_),
56 sg_p_index_map(sg_p_index_map_),
57 sg_g_index_map(sg_g_index_map_),
58 base_g_maps(base_g_maps_),
59 num_p(0),
60 num_g(0),
61 num_p_sg(sg_p_index_map.size()),
62 num_g_sg(sg_g_index_map.size())
63{
64 InArgs me_inargs = me->createInArgs();
65 OutArgs me_outargs = me->createOutArgs();
66 num_p = me_inargs.Np() - num_p_sg;
67 num_g = me_outargs.Ng();
68
69 TEUCHOS_TEST_FOR_EXCEPTION(
70 base_g_maps.size() != num_g_sg, std::logic_error,
71 std::endl
72 << "Error! Stokhos::SGInverseModelEvaluator::SGInverseModelEvaluator():"
73 << " Base response map array size does not match size of index array!");
74}
75
76// Overridden from EpetraExt::ModelEvaluator
77
78Teuchos::RCP<const Epetra_Map>
80get_x_map() const
81{
82 return Teuchos::null;
83}
84
85Teuchos::RCP<const Epetra_Map>
87get_f_map() const
88{
89 return Teuchos::null;
90}
91
92Teuchos::RCP<const Epetra_Map>
94get_p_map(int l) const
95{
96 TEUCHOS_TEST_FOR_EXCEPTION(
97 l >= num_p || l < 0, std::logic_error,
98 std::endl << "Error! Stokhos::SGInverseModelEvaluator::get_p_map():"
99 << " Invalid parameter index l = " << l << std::endl);
100 return me->get_p_map(l);
101}
102
103Teuchos::RCP<const Epetra_Map>
105get_g_map(int l) const
106{
107 TEUCHOS_TEST_FOR_EXCEPTION(
108 l >= num_g || l < 0, std::logic_error,
109 std::endl << "Error! Stokhos::SGInverseModelEvaluator::get_g_map():"
110 << " Invalid response index l = " << l << std::endl);
111 return base_g_maps[l];
112}
113
114Teuchos::RCP<const Teuchos::Array<std::string> >
116get_p_names(int l) const
117{
118 TEUCHOS_TEST_FOR_EXCEPTION(
119 l >= num_p || l < 0, std::logic_error,
120 std::endl << "Error! Stokhos::SGInverseModelEvaluator::get_p_names():"
121 << " Invalid parameter index l = " << l << std::endl);
122 return me->get_p_names(l);
123}
124
125Teuchos::RCP<const Epetra_Vector>
127get_p_init(int l) const
128{
129 TEUCHOS_TEST_FOR_EXCEPTION(
130 l >= num_p || l < 0, std::logic_error,
131 std::endl << "Error! Stokhos::SGInverseModelEvaluator::get_p_init():"
132 << " Invalid parameter index l = " << l << std::endl);
133 return me->get_p_init(l);
134}
135
136EpetraExt::ModelEvaluator::InArgs
138{
139 InArgsSetup inArgs;
140 InArgs me_inargs = me->createInArgs();
141
142 inArgs.setModelEvalDescription(this->description());
143 inArgs.set_Np(num_p);
144 for (int i=0; i<num_p_sg; i++)
145 inArgs.setSupports(IN_ARG_p_sg, sg_p_index_map[i], true);
146 inArgs.setSupports(IN_ARG_sg_basis, me_inargs.supports(IN_ARG_sg_basis));
147 inArgs.setSupports(IN_ARG_sg_quadrature,
148 me_inargs.supports(IN_ARG_sg_quadrature));
149 inArgs.setSupports(IN_ARG_sg_expansion,
150 me_inargs.supports(IN_ARG_sg_expansion));
151
152 return inArgs;
153}
154
155EpetraExt::ModelEvaluator::OutArgs
157{
158 OutArgsSetup outArgs;
159 OutArgs me_outargs = me->createOutArgs();
160
161 outArgs.setModelEvalDescription(this->description());
162
163 outArgs.set_Np_Ng(num_p, num_g);
164 for (int i=0; i<num_g_sg; i++) {
165 outArgs.setSupports(OUT_ARG_g_sg, sg_g_index_map[i], true);
166 for (int j=0; j<num_p; j++)
167 outArgs.setSupports(OUT_ARG_DgDp_sg, sg_g_index_map[i], j,
168 me_outargs.supports(OUT_ARG_DgDp, i, j));
169 }
170
171 return outArgs;
172}
173
174void
176 const OutArgs& outArgs) const
177{
178 // Create underlying inargs
179 InArgs me_inargs = me->createInArgs();
180
181 // Pass parameters
182 for (int i=0; i<num_p; i++)
183 me_inargs.set_p(i, inArgs.get_p(i));
184
185 // Pass SG parameters
186 for (int i=0; i<num_p_sg; i++) {
187 InArgs::sg_const_vector_t p_sg = inArgs.get_p_sg(sg_p_index_map[i]);
188 if (p_sg != Teuchos::null) {
189 me_inargs.set_p(i+num_p, p_sg->getBlockVector());
190 }
191 }
192
193 // Pass SG data
194 if (me_inargs.supports(IN_ARG_sg_basis))
195 me_inargs.set_sg_basis(inArgs.get_sg_basis());
196 if (me_inargs.supports(IN_ARG_sg_quadrature))
197 me_inargs.set_sg_quadrature(inArgs.get_sg_quadrature());
198 if (me_inargs.supports(IN_ARG_sg_expansion))
199 me_inargs.set_sg_expansion(inArgs.get_sg_expansion());
200
201 // Create underlying outargs
202 OutArgs me_outargs = me->createOutArgs();
203
204 // SG Responses
205 for (int i=0; i<num_g_sg; i++) {
206 int ii = sg_g_index_map[i];
207
208 // g
209 OutArgs::sg_vector_t g_sg = outArgs.get_g_sg(ii);
210 if (g_sg != Teuchos::null) {
211 me_outargs.set_g(i, Teuchos::rcp_dynamic_cast<Epetra_Vector>(g_sg->getBlockVector()));
212 }
213
214 // dg/dp
215 for (int j=0; j<num_p; j++) {
216 if (!outArgs.supports(OUT_ARG_DgDp_sg, ii, j).none()) {
217 SGDerivative dgdp_sg = outArgs.get_DgDp_sg(ii,j);
218 Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > dgdp_sg_mv =
219 dgdp_sg.getMultiVector();
220 Teuchos::RCP<Epetra_Operator> dgdp_sg_op =
221 dgdp_sg.getLinearOp();
222 if (dgdp_sg_mv != Teuchos::null) {
223 me_outargs.set_DgDp(
224 i, j, Derivative(dgdp_sg_mv->getBlockMultiVector(),
225 dgdp_sg.getMultiVectorOrientation()));
226 }
227 else if (dgdp_sg_op != Teuchos::null) {
228 me_outargs.set_DgDp(i, j, Derivative(dgdp_sg_op));
229 }
230 }
231 }
232
233 }
234
235 // Compute the functions
236 me->evalModel(me_inargs, me_outargs);
237
238}
Teuchos::RCP< const Epetra_Map > get_g_map(int l) const
Return response map.
Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
Return initial parameters.
int num_g_sg
Number of stochastic response vectors.
Teuchos::RCP< const Epetra_Map > get_x_map() const
Return solution vector map.
SGInverseModelEvaluator(const Teuchos::RCP< EpetraExt::ModelEvaluator > &me, const Teuchos::Array< int > &sg_p_index_map, const Teuchos::Array< int > &sg_g_index_map, const Teuchos::Array< Teuchos::RCP< const Epetra_Map > > &base_g_maps)
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
Return parameter vector map.
Teuchos::RCP< const Epetra_Map > get_f_map() const
Return residual vector map.
int num_p_sg
Number of stochastic parameter vectors.
Teuchos::RCP< EpetraExt::ModelEvaluator > me
Underlying model evaluator.
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
Return array of parameter names.
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > base_g_maps
Base maps of block g vectors.
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Evaluate model on InArgs.