Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_ExplicitModelEvaluator_impl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Panzer: A partial differential equation assembly
4// engine for strongly coupled complex multiphysics systems
5//
6// Copyright 2011 NTESS and the Panzer contributors.
7// SPDX-License-Identifier: BSD-3-Clause
8// *****************************************************************************
9// @HEADER
10
11#ifndef PANZER_EXPLICIT_MODEL_EVALUATOR_IMPL_HPP
12#define PANZER_EXPLICIT_MODEL_EVALUATOR_IMPL_HPP
13
14#include "Thyra_DefaultDiagonalLinearOp.hpp"
15#include "Thyra_LinearOpWithSolveFactoryBase.hpp"
16
17#include "PanzerDiscFE_config.hpp"
18#ifdef PANZER_HAVE_EPETRA_STACK
19#include "Thyra_EpetraModelEvaluator.hpp"
20#include "Thyra_get_Epetra_Operator.hpp"
22#endif
23
24namespace panzer {
25
26template<typename Scalar>
28ExplicitModelEvaluator(const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> > & model,
29 bool constantMassMatrix,
30 bool useLumpedMass,
31 bool applyMassInverse)
32 : Thyra::ModelEvaluatorDelegatorBase<Scalar>(model)
33 , constantMassMatrix_(constantMassMatrix)
34 , massLumping_(useLumpedMass)
35{
36 using Teuchos::RCP;
37 using Teuchos::rcp_dynamic_cast;
38
39 this->applyMassInverse_ = applyMassInverse;
40
41 // extract a panzer::ModelEvaluator if appropriate
42 panzerModel_ = rcp_dynamic_cast<panzer::ModelEvaluator<Scalar> >(model);
43
44#ifdef PANZER_HAVE_EPETRA_STACK
45 // extract a panzer::ModelEvaluator_Epetra if appropriate
46 RCP<Thyra::EpetraModelEvaluator> epME = rcp_dynamic_cast<Thyra::EpetraModelEvaluator>(model);
47 if(epME!=Teuchos::null)
48 panzerEpetraModel_ = rcp_dynamic_cast<const panzer::ModelEvaluator_Epetra>(epME->getEpetraModel());
49#endif
50 // note at this point its possible that panzerModel_ = panzerEpetraModel_ = Teuchos::null
51
53}
54
55template<typename Scalar>
56Thyra::ModelEvaluatorBase::InArgs<Scalar> ExplicitModelEvaluator<Scalar>::
57getNominalValues() const
58{
59 typedef Thyra::ModelEvaluatorBase MEB;
60
61 MEB::InArgs<Scalar> nomVals = createInArgs();
62 nomVals.setArgs(this->getUnderlyingModel()->getNominalValues(),true);
63
64 return nomVals;
65}
66
67template<typename Scalar>
68Thyra::ModelEvaluatorBase::InArgs<Scalar> ExplicitModelEvaluator<Scalar>::
69createInArgs() const
70{
71 return prototypeInArgs_;
72}
73
74template<typename Scalar>
75Thyra::ModelEvaluatorBase::OutArgs<Scalar> ExplicitModelEvaluator<Scalar>::
76createOutArgs() const
77{
78 return prototypeOutArgs_;
79}
80
81template<typename Scalar>
82Teuchos::RCP<panzer::ModelEvaluator<Scalar> > ExplicitModelEvaluator<Scalar>::
84{
85 return Teuchos::rcp_dynamic_cast<panzer::ModelEvaluator<Scalar> >(this->getNonconstUnderlyingModel());
86}
87
88template<typename Scalar>
90evalModelImpl(const Thyra::ModelEvaluatorBase::InArgs<Scalar> &inArgs,
91 const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs) const
92{
93 typedef Thyra::ModelEvaluatorBase MEB;
94 using Teuchos::RCP;
95 RCP<const Thyra::ModelEvaluator<Scalar> > under_me = this->getUnderlyingModel();
96
97 MEB::InArgs<Scalar> under_inArgs = under_me->createInArgs();
98 under_inArgs.setArgs(inArgs);
99
100 // read in the supplied time derivative in case it is needed by the explicit model (e.g. for stabilization) and make sure alpha is set to zero
101 under_inArgs.set_x_dot(inArgs.get_x_dot());
102 under_inArgs.set_alpha(0.0);
103
104 Teuchos::RCP<Thyra::VectorBase<Scalar> > f = outArgs.get_f();
105
106 MEB::OutArgs<Scalar> under_outArgs = under_me->createOutArgs();
107 under_outArgs.setArgs(outArgs);
108 if(f!=Teuchos::null) {
109 // build a scrap vector that will contain the weak residual
110 if(scrap_f_==Teuchos::null)
111 scrap_f_ = Thyra::createMember(*under_me->get_f_space());
112
113 Thyra::assign(scrap_f_.ptr(),0.0);
114 under_outArgs.set_f(scrap_f_);
115 }
116
117 under_me->evalModel(under_inArgs,under_outArgs);
118
119 // build the mass matrix
120 if(invMassMatrix_==Teuchos::null || constantMassMatrix_==false)
121 buildInverseMassMatrix(inArgs);
122
123 if(f!=Teuchos::null)
124 Thyra::Vt_S(scrap_f_.ptr(),-1.0);
125
126 // invert the mass matrix
127 if(f!=Teuchos::null && this->applyMassInverse_) {
128 Thyra::apply(*invMassMatrix_,Thyra::NOTRANS,*scrap_f_,f.ptr());
129 }
130 else if(f!=Teuchos::null){
131 Thyra::V_V(f.ptr(),*scrap_f_);
132 }
133}
134
135template<typename Scalar>
137buildInverseMassMatrix(const Thyra::ModelEvaluatorBase::InArgs<Scalar> &inArgs) const
138{
139 typedef Thyra::ModelEvaluatorBase MEB;
140 using Teuchos::RCP;
141 using Thyra::createMember;
142
143 RCP<const Thyra::ModelEvaluator<Scalar> > me = this->getUnderlyingModel();
144
145 // first allocate space for the mass matrix
146 mass_ = me->create_W_op();
147
148 // intialize a zero to get rid of the x-dot
149 if(zero_==Teuchos::null) {
150 zero_ = Thyra::createMember(*me->get_x_space());
151 Thyra::assign(zero_.ptr(),0.0);
152 }
153
154 // request only the mass matrix from the physics
155 // Model evaluator builds: alpha*u_dot + beta*F(u) = 0
156 MEB::InArgs<Scalar> inArgs_new = me->createInArgs();
157 inArgs_new.setArgs(inArgs);
158 inArgs_new.set_x_dot(inArgs.get_x_dot());
159 inArgs_new.set_alpha(1.0);
160 inArgs_new.set_beta(0.0);
161
162 // set the one time beta to ensure dirichlet conditions
163 // are correctly included in the mass matrix: do it for
164 // both epetra and Tpetra.
165 if(panzerModel_!=Teuchos::null)
166 panzerModel_->setOneTimeDirichletBeta(1.0);
167#ifdef PANZER_HAVE_EPETRA_STACK
168 else if(panzerEpetraModel_!=Teuchos::null)
169 panzerEpetraModel_->setOneTimeDirichletBeta(1.0);
170#endif
171 else {
172 // assuming the underlying model is a delegator, walk through
173 // the decerator hierarchy until you find a panzer::ME or panzer::EpetraME.
174 // If you don't find one, then throw because you are in a load of trouble anyway!
175 setOneTimeDirichletBeta(1.0,*this->getUnderlyingModel());
176 }
177
178 // set only the mass matrix
179 MEB::OutArgs<Scalar> outArgs = me->createOutArgs();
180 outArgs.set_W_op(mass_);
181
182 // this will fill the mass matrix operator
183 me->evalModel(inArgs_new,outArgs);
184
185 // Teuchos::RCP<const Epetra_CrsMatrix> crsMat = Teuchos::rcp_dynamic_cast<const Epetra_CrsMatrix>(Thyra::get_Epetra_Operator(*mass));
186 // EpetraExt::RowMatrixToMatrixMarketFile("expmat.mm",*crsMat);
187
188 if(!massLumping_) {
189 invMassMatrix_ = Thyra::inverse<Scalar>(*me->get_W_factory(),mass_);
190 }
191 else {
192 // build lumped mass matrix (assumes all positive mass entries, does a simple sum)
193 Teuchos::RCP<Thyra::VectorBase<Scalar> > ones = Thyra::createMember(*mass_->domain());
194 Thyra::assign(ones.ptr(),1.0);
195
196 RCP<Thyra::VectorBase<Scalar> > invLumpMass = Thyra::createMember(*mass_->range());
197 Thyra::apply(*mass_,Thyra::NOTRANS,*ones,invLumpMass.ptr());
198 Thyra::reciprocal(*invLumpMass,invLumpMass.ptr());
199
200 invMassMatrix_ = Thyra::diagonal(invLumpMass);
201 }
202}
203
204template<typename Scalar>
207{
208 typedef Thyra::ModelEvaluatorBase MEB;
209
210 MEB::InArgsSetup<Scalar> inArgs(this->getUnderlyingModel()->createInArgs());
211 inArgs.setModelEvalDescription(this->description());
212 inArgs.setSupports(MEB::IN_ARG_alpha,true);
213 inArgs.setSupports(MEB::IN_ARG_beta,true);
214 inArgs.setSupports(MEB::IN_ARG_x_dot,true);
215 prototypeInArgs_ = inArgs;
216
217 MEB::OutArgsSetup<Scalar> outArgs(this->getUnderlyingModel()->createOutArgs());
218 outArgs.setModelEvalDescription(this->description());
219 outArgs.setSupports(MEB::OUT_ARG_W,false);
220 outArgs.setSupports(MEB::OUT_ARG_W_op,false);
221 prototypeOutArgs_ = outArgs;
222}
223
224template<typename Scalar>
226setOneTimeDirichletBeta(double beta,const Thyra::ModelEvaluator<Scalar> & me) const
227{
228 using Teuchos::Ptr;
229 using Teuchos::ptrFromRef;
230 using Teuchos::ptr_dynamic_cast;
231
232 // try to extract base classes that support setOneTimeDirichletBeta
233 Ptr<const panzer::ModelEvaluator<Scalar> > panzerModel = ptr_dynamic_cast<const panzer::ModelEvaluator<Scalar> >(ptrFromRef(me));
234 if(panzerModel!=Teuchos::null) {
235 panzerModel->setOneTimeDirichletBeta(beta);
236 return;
237 }
238 else {
239#ifdef PANZER_HAVE_EPETRA_STACK
240 Ptr<const Thyra::EpetraModelEvaluator> epModel = ptr_dynamic_cast<const Thyra::EpetraModelEvaluator>(ptrFromRef(me));
241 if(epModel!=Teuchos::null) {
242 Ptr<const panzer::ModelEvaluator_Epetra> panzerEpetraModel
243 = ptr_dynamic_cast<const panzer::ModelEvaluator_Epetra>(epModel->getEpetraModel().ptr());
244
245 if(panzerEpetraModel!=Teuchos::null) {
246 panzerEpetraModel->setOneTimeDirichletBeta(beta);
247 return;
248 }
249 }
250#endif
251 }
252
253 // if you get here then the ME is not a panzer::ME or panzer::EpetraME, check
254 // to see if its a delegator
255
256 Ptr<const Thyra::ModelEvaluatorDelegatorBase<Scalar> > delegator
257 = ptr_dynamic_cast<const Thyra::ModelEvaluatorDelegatorBase<Scalar> >(ptrFromRef(me));
258 if(delegator!=Teuchos::null) {
259 setOneTimeDirichletBeta(beta,*delegator->getUnderlyingModel());
260 return;
261 }
262
263 TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
264 "panzer::ExplicitModelEvaluator::setOneTimeDirichletBeta can't find a panzer::ME or a panzer::EpetraME. "
265 "The deepest model is also not a delegator. Thus the recursion failed and an exception was generated.");
266}
267
268} // end namespace panzer
269
270#endif
Teuchos::RCP< const panzer::ModelEvaluator< Scalar > > panzerModel_
Access to the panzer model evaluator pointer (thyra version)
void evalModelImpl(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
Evaluation of model, applies the mass matrix to the RHS.
void setOneTimeDirichletBeta(double beta, const Thyra::ModelEvaluator< Scalar > &me) const
Thyra::ModelEvaluatorBase::InArgs< Scalar > getNominalValues() const
Build the nominal values, modifies the underlying models in args slightly.
void buildInverseMassMatrix(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs) const
Thyra::ModelEvaluatorBase::OutArgs< Scalar > createOutArgs() const
Build the out args, modifies the underlying models in args slightly.
Teuchos::RCP< panzer::ModelEvaluator< Scalar > > getPanzerUnderlyingModel()
Get the underlying panzer::ModelEvaluator.
void buildArgsPrototypes()
Build prototype in/out args.
Thyra::ModelEvaluatorBase::InArgs< Scalar > createInArgs() const
Build the in args, modifies the underlying models in args slightly.
bool applyMassInverse_
Apply mass matrix inverse within the evaluator.