EpetraExt Package Browser (Single Doxygen Collection)
Development
Toggle main menu visibility
Loading...
Searching...
No Matches
example
model_evaluator
DiagonalQuadraticResponseOnlyOpt
EpetraExt_DiagonalQuadraticResponseOnlyModelEvaluator.cpp
Go to the documentation of this file.
1
/*
2
//@HEADER
3
// ***********************************************************************
4
//
5
// EpetraExt: Epetra Extended - Linear Algebra Services Package
6
// Copyright (2011) Sandia Corporation
7
//
8
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9
// the U.S. Government retains certain rights in this software.
10
//
11
// Redistribution and use in source and binary forms, with or without
12
// modification, are permitted provided that the following conditions are
13
// met:
14
//
15
// 1. Redistributions of source code must retain the above copyright
16
// notice, this list of conditions and the following disclaimer.
17
//
18
// 2. Redistributions in binary form must reproduce the above copyright
19
// notice, this list of conditions and the following disclaimer in the
20
// documentation and/or other materials provided with the distribution.
21
//
22
// 3. Neither the name of the Corporation nor the names of the
23
// contributors may be used to endorse or promote products derived from
24
// this software without specific prior written permission.
25
//
26
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37
//
38
// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
39
//
40
// ***********************************************************************
41
//@HEADER
42
*/
43
44
#include "
EpetraExt_DiagonalQuadraticResponseOnlyModelEvaluator.hpp
"
45
#include "Teuchos_ScalarTraits.hpp"
46
#include "
Epetra_SerialComm.h
"
47
#include "
Epetra_CrsMatrix.h
"
48
49
50
namespace
EpetraExt
{
51
52
53
DiagonalQuadraticResponseOnlyModelEvaluator
54
::DiagonalQuadraticResponseOnlyModelEvaluator(
55
const
Teuchos::RCP<Epetra_Comm> &comm,
56
const
int
localDim,
const
double
&pt,
const
double
&p0,
const
double
&scale
57
)
58
:
epetra_comm_
(comm),
scale_
(scale)
59
{
60
61
using
Teuchos::rcp;
62
63
const
int
ng = 1;
64
65
map_p_
= rcp(
new
Epetra_Map
(-1, localDim, 0, *
epetra_comm_
));
66
map_g_
= rcp(
new
Epetra_Map
(ng, ng, 0, *
epetra_comm_
));
67
68
pt_
= rcp(
new
Epetra_Vector
(*
map_p_
));
69
pt_
->PutScalar(pt);
70
71
p0_
= rcp(
new
Epetra_Vector
(*
map_p_
));
72
p0_
->PutScalar(p0);
73
74
}
75
76
77
// Overridden from EpetraExt::ModelEvaluator
78
79
80
Teuchos::RCP<const Epetra_Map>
81
DiagonalQuadraticResponseOnlyModelEvaluator::get_x_map
()
const
82
{
83
return
Teuchos::null;
84
}
85
86
87
Teuchos::RCP<const Epetra_Map>
88
DiagonalQuadraticResponseOnlyModelEvaluator::get_f_map
()
const
89
{
90
return
Teuchos::null;
91
}
92
93
94
Teuchos::RCP<const Epetra_Map>
95
DiagonalQuadraticResponseOnlyModelEvaluator::get_p_map
(
int
l)
const
96
{
97
TEUCHOS_TEST_FOR_EXCEPT(l!=0);
98
return
map_p_
;
99
}
100
101
102
Teuchos::RCP<const Epetra_Map>
103
DiagonalQuadraticResponseOnlyModelEvaluator::get_g_map
(
int
j)
const
104
{
105
TEUCHOS_TEST_FOR_EXCEPT(j!=0);
106
return
map_g_
;
107
}
108
109
110
Teuchos::RCP<const Epetra_Vector>
111
DiagonalQuadraticResponseOnlyModelEvaluator::get_p_init
(
int
l)
const
112
{
113
TEUCHOS_TEST_FOR_EXCEPT(l!=0);
114
return
p0_
;
115
}
116
117
118
EpetraExt::ModelEvaluator::InArgs
119
DiagonalQuadraticResponseOnlyModelEvaluator::createInArgs
()
const
120
{
121
InArgsSetup
inArgs;
122
inArgs.
setModelEvalDescription
(this->description());
123
inArgs.
set_Np
(1);
124
return
inArgs;
125
}
126
127
128
EpetraExt::ModelEvaluator::OutArgs
129
DiagonalQuadraticResponseOnlyModelEvaluator::createOutArgs
()
const
130
{
131
OutArgsSetup
outArgs;
132
outArgs.
setModelEvalDescription
(this->description());
133
outArgs.
set_Np_Ng
(1, 1);
134
outArgs.
setSupports
(
OUT_ARG_DgDp
, 0, 0,
DERIV_TRANS_MV_BY_ROW
);
135
outArgs.
set_DgDp_properties
(
136
0, 0,
DerivativeProperties
(
137
DERIV_LINEARITY_NONCONST
,
138
DERIV_RANK_DEFICIENT
,
139
true
// supportsAdjoint
140
)
141
);
142
return
outArgs;
143
}
144
145
146
void
DiagonalQuadraticResponseOnlyModelEvaluator::evalModel
(
147
const
InArgs
& inArgs,
const
OutArgs
& outArgs
148
)
const
149
{
150
151
using
Teuchos::RCP;
152
using
Teuchos::dyn_cast;
153
using
Teuchos::rcp_dynamic_cast;
154
155
//
156
// Get the input arguments
157
//
158
159
const
Epetra_Vector
&p = *inArgs.
get_p
(0);
160
161
//
162
// Get the output arguments
163
//
164
165
const
RCP<Epetra_Vector> g_out = outArgs.
get_g
(0);
166
167
const
RCP<Epetra_MultiVector> DgDp_trans_out =
168
get_DgDp_mv
(0, 0,outArgs,
DERIV_TRANS_MV_BY_ROW
);
169
170
//
171
// Compute the functions
172
//
173
174
if
(nonnull(g_out) || nonnull(DgDp_trans_out)) {
175
176
Epetra_Vector
p_minus_pt(*
map_p_
);
177
178
p_minus_pt = p;
179
p_minus_pt.
Update
(-1.0, *
pt_
, 1.0);
180
181
if
(nonnull(g_out)) {
182
double
dot[1];
183
p_minus_pt.
Dot
(p_minus_pt, dot);
184
(*g_out)[0] =
scale_
* 0.5 * dot[0];
185
}
186
187
if
(nonnull(DgDp_trans_out)) {
188
(*DgDp_trans_out) = p_minus_pt;
189
DgDp_trans_out->
Scale
(
scale_
);
190
}
191
192
}
193
194
}
195
196
197
}
// namespace EpetraExt
EpetraExt_DiagonalQuadraticResponseOnlyModelEvaluator.hpp
Epetra_CrsMatrix.h
Epetra_SerialComm.h
EpetraExt::DiagonalQuadraticResponseOnlyModelEvaluator::evalModel
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Definition
EpetraExt_DiagonalQuadraticResponseOnlyModelEvaluator.cpp:146
EpetraExt::DiagonalQuadraticResponseOnlyModelEvaluator::createInArgs
InArgs createInArgs() const
Definition
EpetraExt_DiagonalQuadraticResponseOnlyModelEvaluator.cpp:119
EpetraExt::DiagonalQuadraticResponseOnlyModelEvaluator::get_p_map
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
\breif .
Definition
EpetraExt_DiagonalQuadraticResponseOnlyModelEvaluator.cpp:95
EpetraExt::DiagonalQuadraticResponseOnlyModelEvaluator::createOutArgs
OutArgs createOutArgs() const
Definition
EpetraExt_DiagonalQuadraticResponseOnlyModelEvaluator.cpp:129
EpetraExt::DiagonalQuadraticResponseOnlyModelEvaluator::epetra_comm_
Teuchos::RCP< const Epetra_Comm > epetra_comm_
Definition
EpetraExt_DiagonalQuadraticResponseOnlyModelEvaluator.hpp:105
EpetraExt::DiagonalQuadraticResponseOnlyModelEvaluator::get_p_init
Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
Definition
EpetraExt_DiagonalQuadraticResponseOnlyModelEvaluator.cpp:111
EpetraExt::DiagonalQuadraticResponseOnlyModelEvaluator::get_g_map
Teuchos::RCP< const Epetra_Map > get_g_map(int j) const
\breif .
Definition
EpetraExt_DiagonalQuadraticResponseOnlyModelEvaluator.cpp:103
EpetraExt::DiagonalQuadraticResponseOnlyModelEvaluator::map_g_
Teuchos::RCP< const Epetra_Map > map_g_
Definition
EpetraExt_DiagonalQuadraticResponseOnlyModelEvaluator.hpp:107
EpetraExt::DiagonalQuadraticResponseOnlyModelEvaluator::scale_
double scale_
Definition
EpetraExt_DiagonalQuadraticResponseOnlyModelEvaluator.hpp:109
EpetraExt::DiagonalQuadraticResponseOnlyModelEvaluator::get_f_map
Teuchos::RCP< const Epetra_Map > get_f_map() const
Definition
EpetraExt_DiagonalQuadraticResponseOnlyModelEvaluator.cpp:88
EpetraExt::DiagonalQuadraticResponseOnlyModelEvaluator::get_x_map
Teuchos::RCP< const Epetra_Map > get_x_map() const
Definition
EpetraExt_DiagonalQuadraticResponseOnlyModelEvaluator.cpp:81
EpetraExt::DiagonalQuadraticResponseOnlyModelEvaluator::p0_
Teuchos::RCP< Epetra_Vector > p0_
Definition
EpetraExt_DiagonalQuadraticResponseOnlyModelEvaluator.hpp:112
EpetraExt::DiagonalQuadraticResponseOnlyModelEvaluator::pt_
Teuchos::RCP< Epetra_Vector > pt_
Definition
EpetraExt_DiagonalQuadraticResponseOnlyModelEvaluator.hpp:111
EpetraExt::DiagonalQuadraticResponseOnlyModelEvaluator::map_p_
Teuchos::RCP< const Epetra_Map > map_p_
Definition
EpetraExt_DiagonalQuadraticResponseOnlyModelEvaluator.hpp:106
EpetraExt::ModelEvaluator::InArgsSetup
Definition
EpetraExt_ModelEvaluator.h:1291
EpetraExt::ModelEvaluator::InArgsSetup::set_Np
void set_Np(int Np)
Definition
EpetraExt_ModelEvaluator.h:2180
EpetraExt::ModelEvaluator::InArgsSetup::setModelEvalDescription
void setModelEvalDescription(const std::string &modelEvalDescription)
Definition
EpetraExt_ModelEvaluator.h:2174
EpetraExt::ModelEvaluator::InArgs
Definition
EpetraExt_ModelEvaluator.h:135
EpetraExt::ModelEvaluator::InArgs::get_p
Teuchos::RCP< const Epetra_Vector > get_p(int l) const
Definition
EpetraExt_ModelEvaluator.h:1582
EpetraExt::ModelEvaluator::OutArgsSetup
Definition
EpetraExt_ModelEvaluator.h:1306
EpetraExt::ModelEvaluator::OutArgsSetup::setModelEvalDescription
void setModelEvalDescription(const std::string &modelEvalDescription)
Definition
EpetraExt_ModelEvaluator.h:2200
EpetraExt::ModelEvaluator::OutArgsSetup::set_DgDp_properties
void set_DgDp_properties(int j, int l, const DerivativeProperties &properties)
Definition
EpetraExt_ModelEvaluator.h:2313
EpetraExt::ModelEvaluator::OutArgsSetup::set_Np_Ng
void set_Np_Ng(int Np, int Ng)
Definition
EpetraExt_ModelEvaluator.h:2206
EpetraExt::ModelEvaluator::OutArgsSetup::setSupports
void setSupports(EOutArgsMembers arg, bool supports=true)
Definition
EpetraExt_ModelEvaluator.h:2210
EpetraExt::ModelEvaluator::OutArgs
Definition
EpetraExt_ModelEvaluator.h:760
EpetraExt::ModelEvaluator::OutArgs::get_g
Evaluation< Epetra_Vector > get_g(int j) const
Get g(j) where 0 <= j && j < this->Ng().
Definition
EpetraExt_ModelEvaluator.h:1732
EpetraExt::ModelEvaluator::DERIV_TRANS_MV_BY_ROW
@ DERIV_TRANS_MV_BY_ROW
Definition
EpetraExt_ModelEvaluator.h:328
EpetraExt::ModelEvaluator::DERIV_RANK_DEFICIENT
@ DERIV_RANK_DEFICIENT
Definition
EpetraExt_ModelEvaluator.h:413
EpetraExt::ModelEvaluator::OUT_ARG_DgDp
@ OUT_ARG_DgDp
Definition
EpetraExt_ModelEvaluator.h:696
EpetraExt::ModelEvaluator::DERIV_LINEARITY_NONCONST
@ DERIV_LINEARITY_NONCONST
Definition
EpetraExt_ModelEvaluator.h:407
Epetra_Map
Epetra_MultiVector::Scale
int Scale(double ScalarValue)
Epetra_MultiVector::Dot
int Dot(const Epetra_MultiVector &A, double *Result) const
Epetra_MultiVector::Update
int Update(double ScalarA, const Epetra_MultiVector &A, double ScalarThis)
Epetra_Vector
EpetraExt
EpetraExt::BlockCrsMatrix: A class for constructing a distributed block matrix.
Definition
EpetraExt_BlockCrsMatrix.cpp:46
EpetraExt::get_DgDp_mv
Teuchos::RCP< Epetra_MultiVector > get_DgDp_mv(const int j, const int l, const ModelEvaluator::OutArgs &outArgs, const ModelEvaluator::EDerivativeMultiVectorOrientation mvOrientation)
Definition
EpetraExt_ModelEvaluator.cpp:1260
EpetraExt::ModelEvaluator::DerivativeProperties
Definition
EpetraExt_ModelEvaluator.h:417
Generated by
1.17.0