Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_ModelEvaluator_impl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Panzer: A partial differential equation assembly
5 // engine for strongly coupled complex multiphysics systems
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 Roger P. Pawlowski (rppawlo@sandia.gov) and
39 // Eric C. Cyr (eccyr@sandia.gov)
40 // ***********************************************************************
41 // @HEADER
42 
43 #ifndef __Panzer_ModelEvaluator_impl_hpp__
44 #define __Panzer_ModelEvaluator_impl_hpp__
45 
46 #include "Teuchos_DefaultComm.hpp"
47 #include "Teuchos_ArrayRCP.hpp"
48 
49 #include "PanzerDiscFE_config.hpp"
50 #include "Panzer_Traits.hpp"
57 #include "Panzer_GlobalData.hpp"
66 
67 #include "Thyra_TpetraThyraWrappers.hpp"
68 #include "Thyra_SpmdVectorBase.hpp"
69 #include "Thyra_DefaultSpmdVector.hpp"
70 #include "Thyra_DefaultSpmdVectorSpace.hpp"
71 #include "Thyra_DefaultMultiVectorProductVectorSpace.hpp"
72 #include "Thyra_DefaultMultiVectorProductVector.hpp"
73 #include "Thyra_MultiVectorStdOps.hpp"
74 #include "Thyra_VectorStdOps.hpp"
75 
76 #include "Tpetra_CrsMatrix.hpp"
77 
78 // Constructors/Initializers/Accessors
79 
80 template<typename Scalar>
85  const std::vector<Teuchos::RCP<Teuchos::Array<std::string> > >& p_names,
86  const std::vector<Teuchos::RCP<Teuchos::Array<double> > >& p_values,
87  const Teuchos::RCP<const Thyra::LinearOpWithSolveFactoryBase<Scalar> > & solverFactory,
88  const Teuchos::RCP<panzer::GlobalData>& global_data,
89  bool build_transient_support,
90  double t_init)
91  : t_init_(t_init)
92  , num_me_parameters_(0)
93  , do_fd_dfdp_(false)
94  , fd_perturb_size_(1e-7)
95  , require_in_args_refresh_(true)
96  , require_out_args_refresh_(true)
97  , responseLibrary_(rLibrary)
98  , global_data_(global_data)
99  , build_transient_support_(build_transient_support)
100  , lof_(lof)
101  , solverFactory_(solverFactory)
102  , oneTimeDirichletBeta_on_(false)
103  , oneTimeDirichletBeta_(0.0)
104  , build_volume_field_managers_(true)
105  , build_bc_field_managers_(true)
106  , active_evaluation_types_(Sacado::mpl::size<panzer::Traits::EvalTypes>::value, true)
107 {
108  using Teuchos::RCP;
109  using Teuchos::rcp;
110  using Teuchos::rcp_dynamic_cast;
111  using Teuchos::tuple;
112  using Thyra::VectorBase;
113  using Thyra::createMember;
114 
115  TEUCHOS_ASSERT(lof_!=Teuchos::null);
116 
118  ae_tm_.buildObjects(builder);
119 
120  //
121  // Build x, f spaces
122  //
123 
124  // dynamic cast to blocked LOF for now
125  RCP<const ThyraObjFactory<Scalar> > tof = rcp_dynamic_cast<const ThyraObjFactory<Scalar> >(lof,true);
126 
127  x_space_ = tof->getThyraDomainSpace();
128  f_space_ = tof->getThyraRangeSpace();
129 
130  //
131  // Setup parameters
132  //
133  for(std::size_t i=0;i<p_names.size();i++)
134  addParameter(*(p_names[i]),*(p_values[i]));
135 
136  // now that the vector spaces are setup we can allocate the nominal values
137  // (i.e. initial conditions)
139 }
140 
141 template<typename Scalar>
144  const Teuchos::RCP<const Thyra::LinearOpWithSolveFactoryBase<Scalar> > & solverFactory,
145  const Teuchos::RCP<panzer::GlobalData>& global_data,
146  bool build_transient_support,double t_init)
147  : t_init_(t_init)
148  , num_me_parameters_(0)
149  , do_fd_dfdp_(false)
150  , fd_perturb_size_(1e-7)
151  , require_in_args_refresh_(true)
152  , require_out_args_refresh_(true)
153  , global_data_(global_data)
154  , build_transient_support_(build_transient_support)
155  , lof_(lof)
156  , solverFactory_(solverFactory)
157  , oneTimeDirichletBeta_on_(false)
158  , oneTimeDirichletBeta_(0.0)
159  , build_volume_field_managers_(true)
160  , build_bc_field_managers_(true)
161  , active_evaluation_types_(Sacado::mpl::size<panzer::Traits::EvalTypes>::value, true)
162 {
163  using Teuchos::RCP;
164  using Teuchos::rcp_dynamic_cast;
165 
166  TEUCHOS_ASSERT(lof_!=Teuchos::null);
167 
168  //
169  // Build x, f spaces
170  //
171 
172  // dynamic cast to blocked LOF for now
173  RCP<const ThyraObjFactory<Scalar> > tof = rcp_dynamic_cast<const ThyraObjFactory<Scalar> >(lof_,true);
174 
175  x_space_ = tof->getThyraDomainSpace();
176  f_space_ = tof->getThyraRangeSpace();
177 
178  // now that the vector spaces are setup we can allocate the nominal values
179  // (i.e. initial conditions)
181 
182  // allocate a response library so that responses can be added, it will be initialized in "setupModel"
184 }
185 
186 template<typename Scalar>
189 {
190  TEUCHOS_ASSERT(false);
191 }
192 
193 // Public functions overridden from ModelEvaulator
194 
195 template<typename Scalar>
198 {
199  return x_space_;
200 }
201 
202 
203 template<typename Scalar>
206 {
207  return f_space_;
208 }
209 
210 template<typename Scalar>
213 {
214  TEUCHOS_TEST_FOR_EXCEPTION(!(i>=0 && i<num_me_parameters_),std::runtime_error,
215  "panzer::ModelEvaluator::get_p_names: Requested parameter index out of range.");
216 
217  if (i < Teuchos::as<int>(parameters_.size()))
218  return parameters_[i]->names;
219  else if (i < Teuchos::as<int>(parameters_.size()+tangent_space_.size())) {
221  int param_index = i-parameters_.size();
222  std::ostringstream ss;
223  ss << "TANGENT VECTOR: " << param_index;
224  names->push_back(ss.str());
225  return names;
226  }
227  else if (build_transient_support_ && i < Teuchos::as<int>(parameters_.size()+2*tangent_space_.size())) {
229  int param_index = i-parameters_.size()-tangent_space_.size();
230  std::ostringstream ss;
231  ss << "TIME DERIVATIVE TANGENT VECTOR: " << param_index;
232  names->push_back(ss.str());
233  return names;
234  }
235 
236  return Teuchos::null;
237 }
238 
239 template<typename Scalar>
242 {
243  TEUCHOS_TEST_FOR_EXCEPTION(!(i>=0 && i<num_me_parameters_),std::runtime_error,
244  "panzer::ModelEvaluator::get_p_space: Requested parameter index out of range.");
245 
246  if (i < Teuchos::as<int>(parameters_.size()))
247  return parameters_[i]->space;
248  else if (i < Teuchos::as<int>(parameters_.size()+tangent_space_.size()))
249  return tangent_space_[i-parameters_.size()];
250  else if (build_transient_support_ && i < Teuchos::as<int>(parameters_.size()+2*tangent_space_.size()))
251  return tangent_space_[i-parameters_.size()-tangent_space_.size()];
252 
253  return Teuchos::null;
254 }
255 
256 template<typename Scalar>
259 {
260  TEUCHOS_TEST_FOR_EXCEPTION(!(i>=0 && i<Teuchos::as<int>(responses_.size())),std::runtime_error,
261  "panzer::ModelEvaluator::get_g_names: Requested response index out of range.");
262 
263  return Teuchos::ArrayView<const std::string>(&(responses_[i]->name),1);
264 }
265 
266 template<typename Scalar>
267 const std::string &
269 {
270  TEUCHOS_ASSERT(i>=0 &&
271  static_cast<typename std::vector<Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> > >::size_type>(i)<responses_.size());
272 
273  return responses_[i]->name;
274 }
275 
276 template<typename Scalar>
279 {
280  TEUCHOS_ASSERT(i>=0 &&
281  static_cast<typename std::vector<Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> > >::size_type>(i)<responses_.size());
282 
283  return responses_[i]->space;
284 }
285 
286 template<typename Scalar>
287 Thyra::ModelEvaluatorBase::InArgs<Scalar>
289 {
290  return getNominalValues();
291 }
292 
293 template<typename Scalar>
294 Thyra::ModelEvaluatorBase::InArgs<Scalar>
296 {
297  using Teuchos::RCP;
298  using Teuchos::rcp_dynamic_cast;
299 
300  if(require_in_args_refresh_) {
301  typedef Thyra::ModelEvaluatorBase MEB;
302 
303  //
304  // Refresh nominal values, this is primarily adding
305  // new parameters.
306  //
307 
308  MEB::InArgsSetup<Scalar> nomInArgs;
309  nomInArgs = nominalValues_;
310  nomInArgs.setSupports(nominalValues_);
311 
312  // setup parameter support
313  nomInArgs.set_Np(num_me_parameters_);
314  for(std::size_t p=0;p<parameters_.size();p++) {
315  // setup nominal in arguments
316  nomInArgs.set_p(p,parameters_[p]->initial_value);
317 
318  // We explicitly do not set nominal values for tangent parameters
319  // as these are parameters that should be hidden from client code
320  }
321 
322  nominalValues_ = nomInArgs;
323  }
324 
325  // refresh no longer required
326  require_in_args_refresh_ = false;
327 
328  return nominalValues_;
329 }
330 
331 template<typename Scalar>
332 void
334 {
335  typedef Thyra::ModelEvaluatorBase MEB;
336 
337  //
338  // Setup nominal values
339  //
340 
341  MEB::InArgsSetup<Scalar> nomInArgs;
342  nomInArgs.setModelEvalDescription(this->description());
343  nomInArgs.setSupports(MEB::IN_ARG_x);
344  Teuchos::RCP<Thyra::VectorBase<Scalar> > x_nom = Thyra::createMember(x_space_);
345  Thyra::assign(x_nom.ptr(),0.0);
346  nomInArgs.set_x(x_nom);
347  if(build_transient_support_) {
348  nomInArgs.setSupports(MEB::IN_ARG_x_dot,true);
349  nomInArgs.setSupports(MEB::IN_ARG_t,true);
350  nomInArgs.setSupports(MEB::IN_ARG_alpha,true);
351  nomInArgs.setSupports(MEB::IN_ARG_beta,true);
352  nomInArgs.setSupports(MEB::IN_ARG_step_size,true);
353  nomInArgs.setSupports(MEB::IN_ARG_stage_number,true);
354 
355  Teuchos::RCP<Thyra::VectorBase<Scalar> > x_dot_nom = Thyra::createMember(x_space_);
356  Thyra::assign(x_dot_nom.ptr(),0.0);
357  nomInArgs.set_x_dot(x_dot_nom);
358  nomInArgs.set_t(t_init_);
359  nomInArgs.set_alpha(0.0); // these have no meaning initially!
360  nomInArgs.set_beta(0.0);
361  //TODO: is this needed?
362  nomInArgs.set_step_size(0.0);
363  nomInArgs.set_stage_number(1.0);
364  }
365 
366  // setup parameter support -- for each scalar parameter we support the parameter itself and tangent vectors for x, xdot
367  nomInArgs.set_Np(num_me_parameters_);
368  std::size_t v_index = 0;
369  for(std::size_t p=0;p<parameters_.size();p++) {
370  nomInArgs.set_p(p,parameters_[p]->initial_value);
371  if (!parameters_[p]->is_distributed) {
372  Teuchos::RCP<Thyra::VectorBase<Scalar> > v_nom_x = Thyra::createMember(*tangent_space_[v_index]);
373  Thyra::assign(v_nom_x.ptr(),0.0);
374  nomInArgs.set_p(v_index+parameters_.size(),v_nom_x);
375  if (build_transient_support_) {
376  Teuchos::RCP<Thyra::VectorBase<Scalar> > v_nom_xdot = Thyra::createMember(*tangent_space_[v_index]);
377  Thyra::assign(v_nom_xdot.ptr(),0.0);
378  nomInArgs.set_p(v_index+parameters_.size()+tangent_space_.size(),v_nom_xdot);
379  }
380  ++v_index;
381  }
382  }
383 
384  nominalValues_ = nomInArgs;
385 }
386 
387 template <typename Scalar>
389 buildVolumeFieldManagers(const bool value)
390 {
391  build_volume_field_managers_ = value;
392 }
393 
394 template <typename Scalar>
396 buildBCFieldManagers(const bool value)
397 {
398  build_bc_field_managers_ = value;
399 }
400 
401 template <typename Scalar>
404  const std::vector<Teuchos::RCP<panzer::PhysicsBlock> >& physicsBlocks,
405  const std::vector<panzer::BC> & bcs,
406  const panzer::EquationSetFactory & eqset_factory,
407  const panzer::BCStrategyFactory& bc_factory,
410  const Teuchos::ParameterList& closure_models,
411  const Teuchos::ParameterList& user_data,
412  bool writeGraph,const std::string & graphPrefix,
413  const Teuchos::ParameterList& me_params)
414 {
415  // First: build residual assembly engine
417  PANZER_FUNC_TIME_MONITOR_DIFF("panzer::ModelEvaluator::setupModel()",setupModel);
418 
419  {
420  // 1. build Field manager builder
422 
424  {
425  PANZER_FUNC_TIME_MONITOR_DIFF("allocate FieldManagerBuilder",allocFMB);
427  fmb->setActiveEvaluationTypes(active_evaluation_types_);
428  }
429  {
430  PANZER_FUNC_TIME_MONITOR_DIFF("fmb->setWorksetContainer()",setupWorksets);
431  fmb->setWorksetContainer(wc);
432  }
433  if (build_volume_field_managers_) {
434  PANZER_FUNC_TIME_MONITOR_DIFF("fmb->setupVolumeFieldManagers()",setupVolumeFieldManagers);
435  fmb->setupVolumeFieldManagers(physicsBlocks,volume_cm_factory,closure_models,*lof_,user_data);
436  }
437  if (build_bc_field_managers_) {
438  PANZER_FUNC_TIME_MONITOR_DIFF("fmb->setupBCFieldManagers()",setupBCFieldManagers);
439  fmb->setupBCFieldManagers(bcs,physicsBlocks,eqset_factory,bc_cm_factory,bc_factory,closure_models,*lof_,user_data);
440  }
441 
442  // Print Phalanx DAGs
443  if (writeGraph){
444  if (build_volume_field_managers_)
445  fmb->writeVolumeGraphvizDependencyFiles(graphPrefix, physicsBlocks);
446  if (build_bc_field_managers_)
447  fmb->writeBCGraphvizDependencyFiles(graphPrefix+"BC_");
448  }
449 
450  {
451  PANZER_FUNC_TIME_MONITOR_DIFF("AssemblyEngine_TemplateBuilder::buildObjects()",AETM_BuildObjects);
452  panzer::AssemblyEngine_TemplateBuilder builder(fmb,lof_);
453  ae_tm_.buildObjects(builder);
454  }
455  }
456 
457  // Second: build the responses
459 
460  {
461  PANZER_FUNC_TIME_MONITOR_DIFF("build response library",buildResponses);
462 
463  responseLibrary_->initialize(wc,lof_->getRangeGlobalIndexer(),lof_);
464 
465  buildResponses(physicsBlocks,eqset_factory,volume_cm_factory,closure_models,user_data,writeGraph,graphPrefix+"Responses_");
466  buildDistroParamDfDp_RL(wc,physicsBlocks,bcs,eqset_factory,bc_factory,volume_cm_factory,closure_models,user_data,writeGraph,graphPrefix+"Response_DfDp_");
467  buildDistroParamDgDp_RL(wc,physicsBlocks,bcs,eqset_factory,bc_factory,volume_cm_factory,closure_models,user_data,writeGraph,graphPrefix+"Response_DgDp_");
468 
469  do_fd_dfdp_ = false;
470  fd_perturb_size_ = 1.0e-7;
471  if (me_params.isParameter("FD Forward Sensitivities"))
472  do_fd_dfdp_ = me_params.get<bool>("FD Forward Sensitivities");
473  if (me_params.isParameter("FD Perturbation Size"))
474  fd_perturb_size_ = me_params.get<double>("FD Perturbation Size");
475  }
476 }
477 
478 template <typename Scalar>
480 setupAssemblyInArgs(const Thyra::ModelEvaluatorBase::InArgs<Scalar> & inArgs,
481  panzer::AssemblyEngineInArgs & ae_inargs) const
482 {
483  using Teuchos::RCP;
484  using Teuchos::rcp;
485  using Teuchos::rcp_dynamic_cast;
486  using Teuchos::rcp_const_cast;
487  typedef Thyra::ModelEvaluatorBase MEB;
488 
489  // if neccessary build a ghosted container
490  if(Teuchos::is_null(ghostedContainer_)) {
491  ghostedContainer_ = lof_->buildGhostedLinearObjContainer();
492  lof_->initializeGhostedContainer(panzer::LinearObjContainer::X |
495  panzer::LinearObjContainer::Mat, *ghostedContainer_);
496  }
497 
498  bool is_transient = false;
499  if (inArgs.supports(MEB::IN_ARG_x_dot ))
500  is_transient = !Teuchos::is_null(inArgs.get_x_dot());
501 
502  if(Teuchos::is_null(xContainer_))
503  xContainer_ = lof_->buildReadOnlyDomainContainer();
504  if(Teuchos::is_null(xdotContainer_) && is_transient)
505  xdotContainer_ = lof_->buildReadOnlyDomainContainer();
506 
507  const RCP<const Thyra::VectorBase<Scalar> > x = inArgs.get_x();
508  RCP<const Thyra::VectorBase<Scalar> > x_dot; // possibly empty, but otherwise uses x_dot
509 
510  // Make sure construction built in transient support
511  TEUCHOS_TEST_FOR_EXCEPTION(is_transient && !build_transient_support_, std::runtime_error,
512  "ModelEvaluator was not built with transient support enabled!");
513 
514  ae_inargs.container_ = lof_->buildLinearObjContainer(); // we use a new global container
515  ae_inargs.ghostedContainer_ = ghostedContainer_; // we can reuse the ghosted container
516  ae_inargs.alpha = 0.0;
517  ae_inargs.beta = 1.0;
518  ae_inargs.evaluate_transient_terms = false;
519  if (build_transient_support_) {
520  x_dot = inArgs.get_x_dot();
521  ae_inargs.alpha = inArgs.get_alpha();
522  ae_inargs.beta = inArgs.get_beta();
523  ae_inargs.time = inArgs.get_t();
524 
525  ae_inargs.step_size= inArgs.get_step_size();
526  ae_inargs.stage_number = inArgs.get_stage_number();
527  ae_inargs.evaluate_transient_terms = true;
528  }
529 
530  // this member is handled in the individual functions
531  ae_inargs.apply_dirichlet_beta = false;
532 
533  // Set input parameters
534  int num_param_vecs = parameters_.size();
535  for (int i=0; i<num_param_vecs; i++) {
536 
537  RCP<const Thyra::VectorBase<Scalar> > paramVec = inArgs.get_p(i);
538  if ( paramVec!=Teuchos::null && !parameters_[i]->is_distributed) {
539  // non distributed parameters
540 
542  rcp_dynamic_cast<const Thyra::SpmdVectorBase<Scalar> >(paramVec,true)->getLocalData(Teuchos::ptrFromRef(p_data));
543 
544  for (unsigned int j=0; j < parameters_[i]->scalar_value.size(); j++) {
545  parameters_[i]->scalar_value[j].baseValue = p_data[j];
546  parameters_[i]->scalar_value[j].family->setRealValueForAllTypes(parameters_[i]->scalar_value[j].baseValue);
547  }
548  }
549  else if ( paramVec!=Teuchos::null && parameters_[i]->is_distributed) {
550  // distributed parameters
551 
552  std::string key = (*parameters_[i]->names)[0];
553  RCP<GlobalEvaluationData> ged = distrParamGlobalEvaluationData_.getDataObject(key);
554 
555  TEUCHOS_ASSERT(ged!=Teuchos::null);
556 
557  // cast to a LOCPair throwing an exception if the cast doesn't work.
558  RCP<LOCPair_GlobalEvaluationData> loc_pair_ged = rcp_dynamic_cast<LOCPair_GlobalEvaluationData>(ged);
560  if(loc_pair_ged!=Teuchos::null) {
561  // cast to a ThyraObjContainer throwing an exception if the cast doesn't work.
562  RCP<ThyraObjContainer<Scalar> > th_ged = rcp_dynamic_cast<ThyraObjContainer<Scalar> >(loc_pair_ged->getGlobalLOC(),true);
563  th_ged->set_x_th(Teuchos::rcp_const_cast<Thyra::VectorBase<Scalar> >(paramVec));
564  }
565  else {
566  TEUCHOS_ASSERT(ro_ged!=Teuchos::null);
567  ro_ged->setOwnedVector(paramVec);
568  }
569  }
570  }
571 
572  ae_inargs.addGlobalEvaluationData(distrParamGlobalEvaluationData_);
573 
574  // here we are building a container, this operation is fast, simply allocating a struct
575  const RCP<panzer::ThyraObjContainer<Scalar> > thGlobalContainer =
576  Teuchos::rcp_dynamic_cast<panzer::ThyraObjContainer<Scalar> >(ae_inargs.container_);
577 
578  TEUCHOS_ASSERT(!Teuchos::is_null(thGlobalContainer));
579 
580  // Ghosted container objects are zeroed out below only if needed for
581  // a particular calculation. This makes it more efficient than
582  // zeroing out all objects in the container here.
583  // const RCP<panzer::ThyraObjContainer<Scalar> > thGhostedContainer =
584  // Teuchos::rcp_dynamic_cast<panzer::ThyraObjContainer<Scalar> >(ae_inargs.ghostedContainer_);
585 
586  // Set the solution vector (currently all targets require solution).
587  // In the future we may move these into the individual cases below.
588  // A very subtle (and fragile) point: A non-null pointer in global
589  // container triggers export operations during fill. Also, the
590  // introduction of the container is forcing us to cast away const on
591  // arguments that should be const. Another reason to redesign
592  // LinearObjContainer layers.
593  thGlobalContainer->set_x_th(Teuchos::rcp_const_cast<Thyra::VectorBase<Scalar> >(x));
594  xContainer_->setOwnedVector(x);
595  ae_inargs.addGlobalEvaluationData("Solution Gather Container - X",xContainer_);
596 
597  if (is_transient) {
598  thGlobalContainer->set_dxdt_th(Teuchos::rcp_const_cast<Thyra::VectorBase<Scalar> >(x_dot));
599  xdotContainer_->setOwnedVector(x_dot);
600  ae_inargs.addGlobalEvaluationData("Solution Gather Container - Xdot",xdotContainer_);
601  }
602 
603  // Add tangent vectors for x and xdot to GlobalEvaluationData, one for each
604  // scalar parameter vector and parameter within that vector.
605  // Note: The keys for the global evaluation data containers for the tangent
606  // vectors are constructed in EquationSet_AddFieldDefaultImpl::
607  // buildAndRegisterGatherAndOrientationEvaluators().
608  int vIndex(0);
609  for (int i(0); i < num_param_vecs; ++i)
610  {
611  using std::string;
613  using Thyra::VectorBase;
615  if (not parameters_[i]->is_distributed)
616  {
617  auto dxdp = rcp_const_cast<VectorBase<Scalar>>
618  (inArgs.get_p(vIndex + num_param_vecs));
619  if (not dxdp.is_null())
620  {
621  // We need to cast away const because the object container requires
622  // non-const vectors.
623  auto dxdpBlock = rcp_dynamic_cast<ProductVectorBase<Scalar>>(dxdp);
624  int numParams(parameters_[i]->scalar_value.size());
625  for (int j(0); j < numParams; ++j)
626  {
627  RCP<ROVGED> dxdpContainer = lof_->buildReadOnlyDomainContainer();
628  dxdpContainer->setOwnedVector(dxdpBlock->getNonconstVectorBlock(j));
629  string name("X TANGENT GATHER CONTAINER: " +
630  (*parameters_[i]->names)[j]);
631  ae_inargs.addGlobalEvaluationData(name, dxdpContainer);
632  } // end loop over the parameters
633  } // end if (not dxdp.is_null())
634  if (build_transient_support_)
635  {
636  // We need to cast away const because the object container requires
637  // non-const vectors.
638  auto dxdotdp = rcp_const_cast<VectorBase<Scalar>>
639  (inArgs.get_p(vIndex + num_param_vecs + tangent_space_.size()));
640  if (not dxdotdp.is_null())
641  {
642  auto dxdotdpBlock =
643  rcp_dynamic_cast<ProductVectorBase<Scalar>>(dxdotdp);
644  int numParams(parameters_[i]->scalar_value.size());
645  for (int j(0); j < numParams; ++j)
646  {
647  RCP<ROVGED> dxdotdpContainer = lof_->buildReadOnlyDomainContainer();
648  dxdotdpContainer->setOwnedVector(
649  dxdotdpBlock->getNonconstVectorBlock(j));
650  string name("DXDT TANGENT GATHER CONTAINER: " +
651  (*parameters_[i]->names)[j]);
652  ae_inargs.addGlobalEvaluationData(name, dxdotdpContainer);
653  } // end loop over the parameters
654  } // end if (not dxdotdp.is_null())
655  } // end if (build_transient_support_)
656  ++vIndex;
657  } // end if (not parameters_[i]->is_distributed)
658  } // end loop over the parameter vectors
659 } // end of setupAssemblyInArgs()
660 
661 // Private functions overridden from ModelEvaulatorDefaultBase
662 
663 
664 template <typename Scalar>
665 Thyra::ModelEvaluatorBase::OutArgs<Scalar>
667 {
668  typedef Thyra::ModelEvaluatorBase MEB;
669 
670  if(require_out_args_refresh_) {
671  MEB::OutArgsSetup<Scalar> outArgs;
672  outArgs.setModelEvalDescription(this->description());
673  outArgs.set_Np_Ng(num_me_parameters_, responses_.size());
674  outArgs.setSupports(MEB::OUT_ARG_f);
675  outArgs.setSupports(MEB::OUT_ARG_W_op);
676 
677  // add in dg/dx (if appropriate)
678  for(std::size_t i=0;i<responses_.size();i++) {
679  typedef panzer::Traits::Jacobian RespEvalT;
680 
681  // check dg/dx and add it in if appropriate
683  = responseLibrary_->getResponse<RespEvalT>(responses_[i]->name);
684  if(respJacBase!=Teuchos::null) {
685  // cast is guranteed to succeed because of check in addResponse
687  = Teuchos::rcp_dynamic_cast<panzer::ResponseMESupportBase<RespEvalT> >(respJacBase);
688 
689  // class must supppot a derivative
690  if(resp->supportsDerivative()) {
691  outArgs.setSupports(MEB::OUT_ARG_DgDx,i,MEB::DerivativeSupport(MEB::DERIV_MV_GRADIENT_FORM));
692 
693 
694  for(std::size_t p=0;p<parameters_.size();p++) {
695  if(parameters_[p]->is_distributed && parameters_[p]->global_indexer!=Teuchos::null)
696  outArgs.setSupports(MEB::OUT_ARG_DgDp,i,p,MEB::DerivativeSupport(MEB::DERIV_MV_GRADIENT_FORM));
697  if(!parameters_[p]->is_distributed)
698  outArgs.setSupports(MEB::OUT_ARG_DgDp,i,p,MEB::DerivativeSupport(MEB::DERIV_MV_JACOBIAN_FORM));
699  }
700  }
701  }
702  }
703 
704  // setup parameter support
705  for(std::size_t p=0;p<parameters_.size();p++) {
706 
707  if(!parameters_[p]->is_distributed)
708  outArgs.setSupports(MEB::OUT_ARG_DfDp,p,MEB::DerivativeSupport(MEB::DERIV_MV_BY_COL));
709  else if(parameters_[p]->is_distributed && parameters_[p]->global_indexer!=Teuchos::null)
710  outArgs.setSupports(MEB::OUT_ARG_DfDp,p,MEB::DerivativeSupport(MEB::DERIV_LINEAR_OP));
711  }
712 
713  prototypeOutArgs_ = outArgs;
714  }
715 
716  // we don't need to build it anymore
717  require_out_args_refresh_ = false;
718 
719  return prototypeOutArgs_;
720 }
721 
722 template <typename Scalar>
725 create_W_op() const
726 {
728  = Teuchos::rcp_dynamic_cast<const ThyraObjFactory<Scalar> >(lof_,true);
729 
730  return tof->getThyraMatrix();
731 }
732 
733 template <typename Scalar>
737 {
738  return solverFactory_;
739 }
740 
741 template <typename Scalar>
744 create_DfDp_op(int p) const
745 {
746  using Teuchos::RCP;
747  using Teuchos::rcp_dynamic_cast;
748 
749  typedef Thyra::ModelEvaluatorBase MEB;
750 
751  // The code below uses prototypeOutArgs_, so we need to make sure it is
752  // initialized before using it. This happens through createOutArgs(),
753  // however it may be allowable to call create_DfDp_op() before
754  // createOutArgs() is called. Thus we do this here if prototypeOutArgs_
755  // needs to be initialized.
756  //
757  // Previously this was handled in the TEUCHOS_ASSERT below through the call
758  // to Np(), however it isn't a good idea to include code in asserts that is
759  // required for proper execution (the asserts may be removed in an optimized
760  // build, for example).
761  if(require_out_args_refresh_) {
762  this->createOutArgs();
763  }
764 
765  TEUCHOS_ASSERT(0<=p && p<Teuchos::as<int>(parameters_.size()));
766 
767  // assert that DfDp is supported
768  const ParameterObject & po = *parameters_[p];
769 
770  if(po.is_distributed && po.global_indexer!=Teuchos::null) {
771  TEUCHOS_ASSERT(prototypeOutArgs_.supports(MEB::OUT_ARG_DfDp,p).supports(MEB::DERIV_LINEAR_OP));
772 
773  // for a distributed parameter, figure it out from the
774  // response library
775  RCP<Response_Residual<Traits::Jacobian> > response_jacobian
776  = rcp_dynamic_cast<Response_Residual<Traits::Jacobian> >(po.dfdp_rl->template getResponse<Traits::Jacobian>("RESIDUAL"));
777 
778  return response_jacobian->allocateJacobian();
779  }
780  else if(!po.is_distributed) {
781  TEUCHOS_ASSERT(prototypeOutArgs_.supports(MEB::OUT_ARG_DfDp,p).supports(MEB::DERIV_MV_BY_COL));
782 
783  // this is a scalar parameter (its easy to create!)
784  return Thyra::createMember(*get_f_space());
785  }
786 
787  // shourld never get here
788  TEUCHOS_ASSERT(false);
789 
790  return Teuchos::null;
791 }
792 
793 template <typename Scalar>
795 addParameter(const std::string & name,const Scalar & initialValue)
796 {
797  Teuchos::Array<std::string> tmp_names;
798  tmp_names.push_back(name);
799 
800  Teuchos::Array<Scalar> tmp_values;
801  tmp_values.push_back(initialValue);
802 
803  return addParameter(tmp_names,tmp_values);
804 }
805 
806 template <typename Scalar>
809  const Teuchos::Array<Scalar> & initialValues)
810 {
811  using Teuchos::RCP;
812  using Teuchos::rcp;
813  using Teuchos::rcp_dynamic_cast;
814  using Teuchos::ptrFromRef;
815 
816  TEUCHOS_ASSERT(names.size()==initialValues.size());
817 
818  int parameter_index = parameters_.size();
819 
820  // Create parameter object
821  RCP<ParameterObject> param = createScalarParameter(names,initialValues);
822  parameters_.push_back(param);
823 
824  // Create vector space for parameter tangent vector
826  Thyra::multiVectorProductVectorSpace(x_space_, param->names->size());
827  tangent_space_.push_back(tan_space);
828 
829  // The number of model evaluator parameters is the number of model parameters (parameters_.size()) plus a tangent
830  // vector for each scalar parameter (tangent_space_.size()) plus a tangent vector for xdot for each scalar parameter.
831  num_me_parameters_ += 2;
832  if (build_transient_support_)
833  ++num_me_parameters_;
834 
835  require_in_args_refresh_ = true;
836  require_out_args_refresh_ = true;
837  this->resetDefaultBase();
838 
839  return parameter_index;
840 }
841 
842 template <typename Scalar>
844 addDistributedParameter(const std::string & key,
845  const Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> > & vs,
847  const Teuchos::RCP<const Thyra::VectorBase<Scalar> > & initial,
849 {
850  distrParamGlobalEvaluationData_.addDataObject(key,ged);
851 
852  int parameter_index = parameters_.size();
853  parameters_.push_back(createDistributedParameter(key,vs,initial,ugi));
854  ++num_me_parameters_;
855 
856  require_in_args_refresh_ = true;
857  require_out_args_refresh_ = true;
858  this->resetDefaultBase();
859 
860  return parameter_index;
861 }
862 
863 template <typename Scalar>
865 addNonParameterGlobalEvaluationData(const std::string & key,
867 {
868  nonParamGlobalEvaluationData_.addDataObject(key,ged);
869 }
870 
871 template <typename Scalar>
873 addFlexibleResponse(const std::string & responseName,
874  const std::vector<WorksetDescriptor> & wkst_desc,
876 {
877  // add a basic response, use x global indexer to define it
878  builder->setDerivativeInformation(lof_);
879 
880  int respIndex = addResponse(responseName,wkst_desc,*builder);
881 
882  // set the builder for this response
883  responses_[respIndex]->builder = builder;
884 
885  return respIndex;
886 }
887 
888 
889 template <typename Scalar>
892  const Teuchos::RCP<Thyra::VectorBase<Scalar> > & f) const
893 {
894  using Teuchos::RCP;
895  using Teuchos::ArrayRCP;
896  using Teuchos::Array;
897  using Teuchos::tuple;
898  using Teuchos::rcp_dynamic_cast;
899 
900  // if neccessary build a ghosted container
901  if(Teuchos::is_null(ghostedContainer_)) {
902  ghostedContainer_ = lof_->buildGhostedLinearObjContainer();
903  lof_->initializeGhostedContainer(panzer::LinearObjContainer::X |
904  panzer::LinearObjContainer::F,*ghostedContainer_);
905  }
906 
908  ae_inargs.container_ = lof_->buildLinearObjContainer(); // we use a new global container
909  ae_inargs.ghostedContainer_ = ghostedContainer_; // we can reuse the ghosted container
910  ae_inargs.alpha = 0.0;
911  ae_inargs.beta = 1.0;
912  //TODO: is this really needed?
913  ae_inargs.step_size = 0.0;
914  ae_inargs.stage_number = 1.0;
915  ae_inargs.evaluate_transient_terms = false;
916  ae_inargs.addGlobalEvaluationData(nonParamGlobalEvaluationData_);
917  ae_inargs.addGlobalEvaluationData(distrParamGlobalEvaluationData_);
918 
919  // this is the tempory target
920  lof_->initializeContainer(panzer::LinearObjContainer::F,*ae_inargs.container_);
921 
922  // here we are building a container, this operation is fast, simply allocating a struct
923  const RCP<panzer::ThyraObjContainer<Scalar> > thGlobalContainer =
924  Teuchos::rcp_dynamic_cast<panzer::ThyraObjContainer<Scalar> >(ae_inargs.container_);
925 
926  TEUCHOS_ASSERT(!Teuchos::is_null(thGlobalContainer));
927 
928  // Ghosted container objects are zeroed out below only if needed for
929  // a particular calculation. This makes it more efficient than
930  // zeroing out all objects in the container here.
931  const RCP<panzer::ThyraObjContainer<Scalar> > thGhostedContainer =
932  Teuchos::rcp_dynamic_cast<panzer::ThyraObjContainer<Scalar> >(ae_inargs.ghostedContainer_);
933  Thyra::assign(thGhostedContainer->get_f_th().ptr(),0.0);
934 
935  // Set the solution vector (currently all targets require solution).
936  // In the future we may move these into the individual cases below.
937  // A very subtle (and fragile) point: A non-null pointer in global
938  // container triggers export operations during fill. Also, the
939  // introduction of the container is forcing us to cast away const on
940  // arguments that should be const. Another reason to redesign
941  // LinearObjContainer layers.
942  thGlobalContainer->set_x_th(x);
943 
944  // evaluate dirichlet boundary conditions
946  = ae_tm_.template getAsObject<panzer::Traits::Residual>()->evaluateOnlyDirichletBCs(ae_inargs);
947 
948  // allocate the result container
949  RCP<panzer::LinearObjContainer> result = lof_->buildLinearObjContainer(); // we use a new global container
950 
951  // stuff the evaluate boundary conditions into the f spot of the counter ... the x is already filled
952  Teuchos::rcp_dynamic_cast<panzer::ThyraObjContainer<Scalar> >(counter)->set_f_th(
953  thGlobalContainer->get_f_th());
954 
955  // stuff the vector that needs applied dirichlet conditions in the the f spot of the result LOC
956  Teuchos::rcp_dynamic_cast<panzer::ThyraObjContainer<Scalar> >(result)->set_f_th(f);
957 
958  // use the linear object factory to apply the result
959  lof_->applyDirichletBCs(*counter,*result);
960 }
961 
962 template <typename Scalar>
964 evalModel_D2gDx2(int respIndex,
965  const Thyra::ModelEvaluatorBase::InArgs<Scalar> & inArgs,
966  const Teuchos::RCP<const Thyra::VectorBase<Scalar> > & delta_x,
967  const Teuchos::RCP<Thyra::VectorBase<Scalar> > & D2gDx2) const
968 {
969 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
970 
971  // set model parameters from supplied inArgs
972  setParameters(inArgs);
973 
974  {
975  std::string responseName = responses_[respIndex]->name;
977  = Teuchos::rcp_dynamic_cast<panzer::ResponseMESupportBase<panzer::Traits::Hessian> >(
978  responseLibrary_->getResponse<panzer::Traits::Hessian>(responseName));
979  resp->setDerivative(D2gDx2);
980  }
981 
982  // setup all the assembly in arguments (this is parameters and
983  // x/x_dot). At this point with the exception of the one time dirichlet
984  // beta that is all thats neccessary.
986  setupAssemblyInArgs(inArgs,ae_inargs);
987 
988  ae_inargs.beta = 1.0;
989 
990  auto deltaXContainer = lof_->buildReadOnlyDomainContainer();
991  deltaXContainer->setOwnedVector(delta_x);
992  ae_inargs.addGlobalEvaluationData("DELTA_Solution Gather Container",deltaXContainer);
993 
994  // evaluate responses
995  responseLibrary_->addResponsesToInArgs<panzer::Traits::Hessian>(ae_inargs);
996  responseLibrary_->evaluate<panzer::Traits::Hessian>(ae_inargs);
997 
998  // reset parameters back to nominal values
999  resetParameters();
1000 #else
1001  (void)respIndex;
1002  (void)inArgs;
1003  (void)delta_x;
1004  (void)D2gDx2;
1005  TEUCHOS_ASSERT(false);
1006 #endif
1007 }
1008 
1009 template <typename Scalar>
1011 evalModel_D2gDxDp(int respIndex,
1012  int pIndex,
1013  const Thyra::ModelEvaluatorBase::InArgs<Scalar> & inArgs,
1014  const Teuchos::RCP<const Thyra::VectorBase<Scalar> > & delta_p,
1015  const Teuchos::RCP<Thyra::VectorBase<Scalar> > & D2gDxDp) const
1016 {
1017 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
1018 
1019  // set model parameters from supplied inArgs
1020  setParameters(inArgs);
1021 
1022  {
1023  std::string responseName = responses_[respIndex]->name;
1025  = Teuchos::rcp_dynamic_cast<panzer::ResponseMESupportBase<panzer::Traits::Hessian> >(
1026  responseLibrary_->getResponse<panzer::Traits::Hessian>(responseName));
1027  resp->setDerivative(D2gDxDp);
1028  }
1029 
1030  // setup all the assembly in arguments (this is parameters and
1031  // x/x_dot). At this point with the exception of the one time dirichlet
1032  // beta that is all thats neccessary.
1033  panzer::AssemblyEngineInArgs ae_inargs;
1034  setupAssemblyInArgs(inArgs,ae_inargs);
1035 
1036  ae_inargs.beta = 1.0;
1037  ae_inargs.second_sensitivities_name = (*parameters_[pIndex]->names)[0]; // distributed parameters have one name!
1038 
1039  auto deltaPContainer = parameters_[pIndex]->dfdp_rl->getLinearObjFactory()->buildReadOnlyDomainContainer();
1040  deltaPContainer->setOwnedVector(delta_p);
1041  ae_inargs.addGlobalEvaluationData("DELTA_"+(*parameters_[pIndex]->names)[0],deltaPContainer);
1042 
1043  // evaluate responses
1044  responseLibrary_->addResponsesToInArgs<panzer::Traits::Hessian>(ae_inargs);
1045  responseLibrary_->evaluate<panzer::Traits::Hessian>(ae_inargs);
1046 
1047  // reset parameters back to nominal values
1048  resetParameters();
1049 #else
1050  (void)respIndex;
1051  (void)pIndex;
1052  (void)inArgs;
1053  (void)delta_p;
1054  (void)D2gDxDp;
1055  TEUCHOS_ASSERT(false);
1056 #endif
1057 }
1058 
1059 template <typename Scalar>
1061 evalModel_D2gDp2(int respIndex,
1062  int pIndex,
1063  const Thyra::ModelEvaluatorBase::InArgs<Scalar> & inArgs,
1064  const Teuchos::RCP<const Thyra::VectorBase<Scalar> > & delta_p,
1065  const Teuchos::RCP<Thyra::VectorBase<Scalar> > & D2gDp2) const
1066 {
1067 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
1068 
1069  // set model parameters from supplied inArgs
1070  setParameters(inArgs);
1071 
1072  ResponseLibrary<Traits> & rLibrary = *parameters_[pIndex]->dgdp_rl;
1073 
1074  {
1075  std::string responseName = responses_[respIndex]->name;
1077  = Teuchos::rcp_dynamic_cast<panzer::ResponseMESupportBase<panzer::Traits::Hessian> >(
1078  rLibrary.getResponse<panzer::Traits::Hessian>(responseName));
1079  resp->setDerivative(D2gDp2);
1080  }
1081 
1082  // setup all the assembly in arguments (this is parameters and
1083  // x/x_dot). At this point with the exception of the one time dirichlet
1084  // beta that is all thats neccessary.
1085  panzer::AssemblyEngineInArgs ae_inargs;
1086  setupAssemblyInArgs(inArgs,ae_inargs);
1087 
1088  ae_inargs.gather_seeds.push_back(1.0); // this assumes that gather point is always the zero index of
1089  // gather seeds
1090  ae_inargs.first_sensitivities_name = (*parameters_[pIndex]->names)[0]; // distributed parameters have one name!
1091  ae_inargs.second_sensitivities_name = (*parameters_[pIndex]->names)[0]; // distributed parameters have one name!
1092 
1093  auto deltaPContainer = parameters_[pIndex]->dfdp_rl->getLinearObjFactory()->buildReadOnlyDomainContainer();
1094  deltaPContainer->setOwnedVector(delta_p);
1095  ae_inargs.addGlobalEvaluationData("DELTA_"+(*parameters_[pIndex]->names)[0],deltaPContainer);
1096 
1097  // evaluate responses
1098  rLibrary.addResponsesToInArgs<panzer::Traits::Hessian>(ae_inargs);
1099  rLibrary.evaluate<panzer::Traits::Hessian>(ae_inargs);
1100 
1101  // reset parameters back to nominal values
1102  resetParameters();
1103 #else
1104  (void)respIndex;
1105  (void)pIndex;
1106  (void)inArgs;
1107  (void)delta_p;
1108  (void)D2gDp2;
1109  TEUCHOS_ASSERT(false);
1110 #endif
1111 }
1112 
1113 template <typename Scalar>
1115 evalModel_D2gDpDx(int respIndex,
1116  int pIndex,
1117  const Thyra::ModelEvaluatorBase::InArgs<Scalar> & inArgs,
1118  const Teuchos::RCP<const Thyra::VectorBase<Scalar> > & delta_x,
1119  const Teuchos::RCP<Thyra::VectorBase<Scalar> > & D2gDpDx) const
1120 {
1121 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
1122 
1123  // set model parameters from supplied inArgs
1124  setParameters(inArgs);
1125 
1126  ResponseLibrary<Traits> & rLibrary = *parameters_[pIndex]->dgdp_rl;
1127 
1128  {
1129  std::string responseName = responses_[respIndex]->name;
1131  = Teuchos::rcp_dynamic_cast<panzer::ResponseMESupportBase<panzer::Traits::Hessian> >(
1132  rLibrary.getResponse<panzer::Traits::Hessian>(responseName));
1133  resp->setDerivative(D2gDpDx);
1134  }
1135 
1136  // setup all the assembly in arguments (this is parameters and
1137  // x/x_dot). At this point with the exception of the one time dirichlet
1138  // beta that is all thats neccessary.
1139  panzer::AssemblyEngineInArgs ae_inargs;
1140  setupAssemblyInArgs(inArgs,ae_inargs);
1141 
1142  ae_inargs.gather_seeds.push_back(1.0); // this assumes that gather point is always the zero index of
1143  // gather seeds
1144  ae_inargs.first_sensitivities_name = (*parameters_[pIndex]->names)[0]; // distributed parameters have one name!
1145  ae_inargs.second_sensitivities_name = "";
1146 
1147  auto deltaXContainer = lof_->buildReadOnlyDomainContainer();
1148  deltaXContainer->setOwnedVector(delta_x);
1149  ae_inargs.addGlobalEvaluationData("DELTA_Solution Gather Container",deltaXContainer);
1150 
1151  // evaluate responses
1152  rLibrary.addResponsesToInArgs<panzer::Traits::Hessian>(ae_inargs);
1153  rLibrary.evaluate<panzer::Traits::Hessian>(ae_inargs);
1154 
1155  // reset parameters back to nominal values
1156  resetParameters();
1157 #else
1158  (void)respIndex;
1159  (void)pIndex;
1160  (void)inArgs;
1161  (void)delta_x;
1162  (void)D2gDpDx;
1163  TEUCHOS_ASSERT(false);
1164 #endif
1165 }
1166 
1167 template <typename Scalar>
1169 evalModel_D2fDx2(const Thyra::ModelEvaluatorBase::InArgs<Scalar> & inArgs,
1170  const Teuchos::RCP<const Thyra::VectorBase<Scalar> > & delta_x,
1171  const Teuchos::RCP<Thyra::LinearOpBase<Scalar> > & D2fDx2) const
1172 {
1173 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
1174 
1175  using Teuchos::RCP;
1176  using Teuchos::ArrayRCP;
1177  using Teuchos::Array;
1178  using Teuchos::tuple;
1179  using Teuchos::rcp_dynamic_cast;
1180 
1181  typedef Thyra::ModelEvaluatorBase MEB;
1182 
1183  // Transient or steady-state evaluation is determined by the x_dot
1184  // vector. If this RCP is null, then we are doing a steady-state
1185  // fill.
1186  bool is_transient = false;
1187  if (inArgs.supports(MEB::IN_ARG_x_dot ))
1188  is_transient = !Teuchos::is_null(inArgs.get_x_dot());
1189 
1190  // Make sure construction built in transient support
1191  TEUCHOS_TEST_FOR_EXCEPTION(is_transient && !build_transient_support_, std::runtime_error,
1192  "ModelEvaluator was not built with transient support enabled!");
1193 
1194  //
1195  // Get the output arguments
1196  //
1197  const RCP<Thyra::LinearOpBase<Scalar> > W_out = D2fDx2;
1198 
1199  // setup all the assembly in arguments (this is parameters and
1200  // x/x_dot). At this point with the exception of the one time dirichlet
1201  // beta that is all thats neccessary.
1202  panzer::AssemblyEngineInArgs ae_inargs;
1203  setupAssemblyInArgs(inArgs,ae_inargs);
1204 
1205  auto deltaXContainer = lof_->buildReadOnlyDomainContainer();
1206  deltaXContainer->setOwnedVector(delta_x);
1207  ae_inargs.addGlobalEvaluationData("DELTA_Solution Gather Container",deltaXContainer);
1208 
1209  // set model parameters from supplied inArgs
1210  setParameters(inArgs);
1211 
1212  // handle application of the one time dirichlet beta in the
1213  // assembly engine. Note that this has to be set explicitly
1214  // each time because this badly breaks encapsulation. Essentially
1215  // we must work around the model evaluator abstraction!
1216  if(oneTimeDirichletBeta_on_) {
1217  ae_inargs.dirichlet_beta = oneTimeDirichletBeta_;
1218  ae_inargs.apply_dirichlet_beta = true;
1219 
1220  oneTimeDirichletBeta_on_ = false;
1221  }
1222 
1223  // here we are building a container, this operation is fast, simply allocating a struct
1224  const RCP<panzer::ThyraObjContainer<Scalar> > thGlobalContainer =
1225  Teuchos::rcp_dynamic_cast<panzer::ThyraObjContainer<Scalar> >(ae_inargs.container_);
1226  const RCP<panzer::ThyraObjContainer<Scalar> > thGhostedContainer =
1227  Teuchos::rcp_dynamic_cast<panzer::ThyraObjContainer<Scalar> >(ae_inargs.ghostedContainer_);
1228 
1229  {
1230  PANZER_FUNC_TIME_MONITOR("panzer::ModelEvaluator::evalModel(D2fDx2)");
1231 
1232  // this dummy nonsense is needed only for scattering dirichlet conditions
1233  RCP<Thyra::VectorBase<Scalar> > dummy_f = Thyra::createMember(f_space_);
1234  thGlobalContainer->set_f_th(dummy_f);
1235  thGlobalContainer->set_A_th(W_out);
1236 
1237  // Zero values in ghosted container objects
1238  thGhostedContainer->initializeMatrix(0.0);
1239 
1240  ae_tm_.template getAsObject<panzer::Traits::Hessian>()->evaluate(ae_inargs);
1241  }
1242 
1243  // HACK: set A to null before calling responses to avoid touching the
1244  // the Jacobian after it has been properly assembled. Should be fixed
1245  // by using a modified version of ae_inargs instead.
1246  thGlobalContainer->set_A_th(Teuchos::null);
1247 
1248  // Holding a rcp to f produces a seg fault in Rythmos when the next
1249  // f comes in and the resulting dtor is called. Need to discuss
1250  // with Ross. Clearing all references here works!
1251 
1252  thGlobalContainer->set_x_th(Teuchos::null);
1253  thGlobalContainer->set_dxdt_th(Teuchos::null);
1254  thGlobalContainer->set_f_th(Teuchos::null);
1255  thGlobalContainer->set_A_th(Teuchos::null);
1256 
1257  // reset parameters back to nominal values
1258  resetParameters();
1259 #else
1260  (void)inArgs;
1261  (void)delta_x;
1262  (void)D2fDx2;
1263  TEUCHOS_ASSERT(false);
1264 #endif
1265 }
1266 
1267 template <typename Scalar>
1270  const Thyra::ModelEvaluatorBase::InArgs<Scalar> & inArgs,
1271  const Teuchos::RCP<const Thyra::VectorBase<Scalar> > & delta_p,
1272  const Teuchos::RCP<Thyra::LinearOpBase<Scalar> > & D2fDxDp) const
1273 {
1274 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
1275 
1276  using Teuchos::RCP;
1277  using Teuchos::ArrayRCP;
1278  using Teuchos::Array;
1279  using Teuchos::tuple;
1280  using Teuchos::rcp_dynamic_cast;
1281 
1282  typedef Thyra::ModelEvaluatorBase MEB;
1283 
1284  // Transient or steady-state evaluation is determined by the x_dot
1285  // vector. If this RCP is null, then we are doing a steady-state
1286  // fill.
1287  bool is_transient = false;
1288  if (inArgs.supports(MEB::IN_ARG_x_dot ))
1289  is_transient = !Teuchos::is_null(inArgs.get_x_dot());
1290 
1291  // Make sure construction built in transient support
1292  TEUCHOS_TEST_FOR_EXCEPTION(is_transient && !build_transient_support_, std::runtime_error,
1293  "ModelEvaluator was not built with transient support enabled!");
1294 
1295  //
1296  // Get the output arguments
1297  //
1298  const RCP<Thyra::LinearOpBase<Scalar> > W_out = D2fDxDp;
1299 
1300  // setup all the assembly in arguments (this is parameters and
1301  // x/x_dot). At this point with the exception of the one time dirichlet
1302  // beta that is all thats neccessary.
1303  panzer::AssemblyEngineInArgs ae_inargs;
1304  setupAssemblyInArgs(inArgs,ae_inargs);
1305 
1306  ae_inargs.second_sensitivities_name = (*parameters_[pIndex]->names)[0]; // distributed parameters have one name!
1307 
1308  auto deltaPContainer = parameters_[pIndex]->dfdp_rl->getLinearObjFactory()->buildReadOnlyDomainContainer();
1309  deltaPContainer->setOwnedVector(delta_p);
1310  ae_inargs.addGlobalEvaluationData("DELTA_"+(*parameters_[pIndex]->names)[0],deltaPContainer);
1311 
1312  // set model parameters from supplied inArgs
1313  setParameters(inArgs);
1314 
1315  // handle application of the one time dirichlet beta in the
1316  // assembly engine. Note that this has to be set explicitly
1317  // each time because this badly breaks encapsulation. Essentially
1318  // we must work around the model evaluator abstraction!
1319  if(oneTimeDirichletBeta_on_) {
1320  ae_inargs.dirichlet_beta = oneTimeDirichletBeta_;
1321  ae_inargs.apply_dirichlet_beta = true;
1322 
1323  oneTimeDirichletBeta_on_ = false;
1324  }
1325 
1326  // here we are building a container, this operation is fast, simply allocating a struct
1327  const RCP<panzer::ThyraObjContainer<Scalar> > thGlobalContainer =
1328  Teuchos::rcp_dynamic_cast<panzer::ThyraObjContainer<Scalar> >(ae_inargs.container_);
1329  const RCP<panzer::ThyraObjContainer<Scalar> > thGhostedContainer =
1330  Teuchos::rcp_dynamic_cast<panzer::ThyraObjContainer<Scalar> >(ae_inargs.ghostedContainer_);
1331 
1332  {
1333  PANZER_FUNC_TIME_MONITOR("panzer::ModelEvaluator::evalModel(D2fDxDp)");
1334 
1335  // this dummy nonsense is needed only for scattering dirichlet conditions
1336  RCP<Thyra::VectorBase<Scalar> > dummy_f = Thyra::createMember(f_space_);
1337  thGlobalContainer->set_f_th(dummy_f);
1338  thGlobalContainer->set_A_th(W_out);
1339 
1340  // Zero values in ghosted container objects
1341  thGhostedContainer->initializeMatrix(0.0);
1342 
1343  ae_tm_.template getAsObject<panzer::Traits::Hessian>()->evaluate(ae_inargs);
1344  }
1345 
1346  // HACK: set A to null before calling responses to avoid touching the
1347  // the Jacobian after it has been properly assembled. Should be fixed
1348  // by using a modified version of ae_inargs instead.
1349  thGlobalContainer->set_A_th(Teuchos::null);
1350 
1351  // Holding a rcp to f produces a seg fault in Rythmos when the next
1352  // f comes in and the resulting dtor is called. Need to discuss
1353  // with Ross. Clearing all references here works!
1354 
1355  thGlobalContainer->set_x_th(Teuchos::null);
1356  thGlobalContainer->set_dxdt_th(Teuchos::null);
1357  thGlobalContainer->set_f_th(Teuchos::null);
1358  thGlobalContainer->set_A_th(Teuchos::null);
1359 
1360  // reset parameters back to nominal values
1361  resetParameters();
1362 #else
1363  (void)pIndex;
1364  (void)inArgs;
1365  (void)delta_p;
1366  (void)D2fDxDp;
1367  TEUCHOS_ASSERT(false);
1368 #endif
1369 }
1370 
1371 template <typename Scalar>
1374  const Thyra::ModelEvaluatorBase::InArgs<Scalar> & inArgs,
1375  const Teuchos::RCP<const Thyra::VectorBase<Scalar> > & delta_x,
1376  const Teuchos::RCP<Thyra::LinearOpBase<Scalar> > & D2fDpDx) const
1377 {
1378 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
1379  using Teuchos::RCP;
1380  using Teuchos::rcp_dynamic_cast;
1381  using Teuchos::null;
1382 
1383  // parameter is not distributed
1384  TEUCHOS_ASSERT(parameters_[pIndex]->is_distributed);
1385 
1386  // parameter is distributed but has no global indexer.
1387  // thus the user doesn't want sensitivities!
1388  TEUCHOS_ASSERT(parameters_[pIndex]->dfdp_rl!=null);
1389 
1390  ResponseLibrary<Traits> & rLibrary = *parameters_[pIndex]->dfdp_rl;
1391 
1392  // get the response and tell it to fill the derivative operator
1393  RCP<Response_Residual<Traits::Hessian> > response_hessian =
1394  rcp_dynamic_cast<Response_Residual<Traits::Hessian> >(rLibrary.getResponse<Traits::Hessian>("RESIDUAL"));
1395  response_hessian->setHessian(D2fDpDx);
1396 
1397  // setup all the assembly in arguments (this is parameters and x/x_dot).
1398  // make sure the correct seeding is performed
1399  panzer::AssemblyEngineInArgs ae_inargs;
1400  setupAssemblyInArgs(inArgs,ae_inargs);
1401 
1402  auto deltaXContainer = lof_->buildReadOnlyDomainContainer();
1403  deltaXContainer->setOwnedVector(delta_x);
1404  ae_inargs.addGlobalEvaluationData("DELTA_Solution Gather Container",deltaXContainer);
1405 
1406  ae_inargs.gather_seeds.push_back(1.0); // this assumes that gather point is always the zero index of
1407  // gather seeds
1408  ae_inargs.first_sensitivities_name = (*parameters_[pIndex]->names)[0]; // distributed parameters have one name!
1409  ae_inargs.second_sensitivities_name = "";
1410 
1411  rLibrary.addResponsesToInArgs<Traits::Hessian>(ae_inargs);
1412  rLibrary.evaluate<Traits::Hessian>(ae_inargs);
1413 #else
1414  (void)pIndex;
1415  (void)inArgs;
1416  (void)delta_x;
1417  (void)D2fDpDx;
1418  TEUCHOS_ASSERT(false);
1419 #endif
1420 }
1421 
1422 template <typename Scalar>
1424 evalModel_D2fDp2(int pIndex,
1425  const Thyra::ModelEvaluatorBase::InArgs<Scalar> & inArgs,
1426  const Teuchos::RCP<const Thyra::VectorBase<Scalar> > & delta_p,
1427  const Teuchos::RCP<Thyra::LinearOpBase<Scalar> > & D2fDp2) const
1428 {
1429 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
1430  using Teuchos::RCP;
1431  using Teuchos::rcp_dynamic_cast;
1432  using Teuchos::null;
1433 
1434  // parameter is not distributed
1435  TEUCHOS_ASSERT(parameters_[pIndex]->is_distributed);
1436 
1437  // parameter is distributed but has no global indexer.
1438  // thus the user doesn't want sensitivities!
1439  TEUCHOS_ASSERT(parameters_[pIndex]->dfdp_rl!=null);
1440 
1441  ResponseLibrary<Traits> & rLibrary = *parameters_[pIndex]->dfdp_rl;
1442 
1443  // get the response and tell it to fill the derivative operator
1444  RCP<Response_Residual<Traits::Hessian> > response_hessian =
1445  rcp_dynamic_cast<Response_Residual<Traits::Hessian> >(rLibrary.getResponse<Traits::Hessian>("RESIDUAL"));
1446  response_hessian->setHessian(D2fDp2);
1447 
1448  // setup all the assembly in arguments (this is parameters and x/x_dot).
1449  // make sure the correct seeding is performed
1450  panzer::AssemblyEngineInArgs ae_inargs;
1451  setupAssemblyInArgs(inArgs,ae_inargs);
1452 
1453  auto deltaPContainer = parameters_[pIndex]->dfdp_rl->getLinearObjFactory()->buildReadOnlyDomainContainer();
1454  deltaPContainer->setOwnedVector(delta_p);
1455  ae_inargs.addGlobalEvaluationData("DELTA_"+(*parameters_[pIndex]->names)[0],deltaPContainer);
1456 
1457  ae_inargs.gather_seeds.push_back(1.0); // this assumes that gather point is always the zero index of
1458  // gather seeds
1459  ae_inargs.first_sensitivities_name = (*parameters_[pIndex]->names)[0]; // distributed parameters have one name!
1460  ae_inargs.second_sensitivities_name = (*parameters_[pIndex]->names)[0]; // distributed parameters have one name!
1461 
1462  rLibrary.addResponsesToInArgs<Traits::Hessian>(ae_inargs);
1463  rLibrary.evaluate<Traits::Hessian>(ae_inargs);
1464 #else
1465  (void)pIndex;
1466  (void)inArgs;
1467  (void)delta_p;
1468  (void)D2fDp2;
1469  TEUCHOS_ASSERT(false);
1470 #endif
1471 }
1472 
1473 template <typename Scalar>
1475 evalModelImpl(const Thyra::ModelEvaluatorBase::InArgs<Scalar> &inArgs,
1476  const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs) const
1477 {
1478  evalModelImpl_basic(inArgs,outArgs);
1479 
1480  // evaluate responses...uses the stored assembly arguments and containers
1481  if(required_basic_g(outArgs))
1482  evalModelImpl_basic_g(inArgs,outArgs);
1483 
1484  // evaluate response derivatives
1485  if(required_basic_dgdx(outArgs))
1486  evalModelImpl_basic_dgdx(inArgs,outArgs);
1487 
1488  // evaluate response derivatives to scalar parameters
1489  if(required_basic_dgdp_scalar(outArgs))
1490  evalModelImpl_basic_dgdp_scalar(inArgs,outArgs);
1491 
1492  // evaluate response derivatives to distributed parameters
1493  if(required_basic_dgdp_distro(outArgs))
1494  evalModelImpl_basic_dgdp_distro(inArgs,outArgs);
1495 
1496  if(required_basic_dfdp_scalar(outArgs)) {
1497  if (do_fd_dfdp_)
1498  evalModelImpl_basic_dfdp_scalar_fd(inArgs,outArgs);
1499  else
1500  evalModelImpl_basic_dfdp_scalar(inArgs,outArgs);
1501  }
1502 
1503  if(required_basic_dfdp_distro(outArgs))
1504  evalModelImpl_basic_dfdp_distro(inArgs,outArgs);
1505 }
1506 
1507 template <typename Scalar>
1509 evalModelImpl_basic(const Thyra::ModelEvaluatorBase::InArgs<Scalar> &inArgs,
1510  const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs) const
1511 {
1512  using Teuchos::RCP;
1513  using Teuchos::ArrayRCP;
1514  using Teuchos::Array;
1515  using Teuchos::tuple;
1516  using Teuchos::rcp_dynamic_cast;
1517 
1518  typedef Thyra::ModelEvaluatorBase MEB;
1519 
1520  // Transient or steady-state evaluation is determined by the x_dot
1521  // vector. If this RCP is null, then we are doing a steady-state
1522  // fill.
1523  bool is_transient = false;
1524  if (inArgs.supports(MEB::IN_ARG_x_dot ))
1525  is_transient = !Teuchos::is_null(inArgs.get_x_dot());
1526 
1527  // Make sure construction built in transient support
1528  TEUCHOS_TEST_FOR_EXCEPTION(is_transient && !build_transient_support_, std::runtime_error,
1529  "ModelEvaluator was not built with transient support enabled!");
1530 
1531  //
1532  // Get the output arguments
1533  //
1534  const RCP<Thyra::VectorBase<Scalar> > f_out = outArgs.get_f();
1535  const RCP<Thyra::LinearOpBase<Scalar> > W_out = outArgs.get_W_op();
1536 
1537  // see if the user wants us to do anything
1538  if(Teuchos::is_null(f_out) && Teuchos::is_null(W_out) ) {
1539  return;
1540  }
1541 
1542  // setup all the assembly in arguments (this is parameters and
1543  // x/x_dot). At this point with the exception of the one time dirichlet
1544  // beta that is all thats neccessary.
1545  panzer::AssemblyEngineInArgs ae_inargs;
1546  setupAssemblyInArgs(inArgs,ae_inargs);
1547 
1548  // set model parameters from supplied inArgs
1549  setParameters(inArgs);
1550 
1551  // handle application of the one time dirichlet beta in the
1552  // assembly engine. Note that this has to be set explicitly
1553  // each time because this badly breaks encapsulation. Essentially
1554  // we must work around the model evaluator abstraction!
1555  if(oneTimeDirichletBeta_on_) {
1556  ae_inargs.dirichlet_beta = oneTimeDirichletBeta_;
1557  ae_inargs.apply_dirichlet_beta = true;
1558 
1559  oneTimeDirichletBeta_on_ = false;
1560  }
1561 
1562  // here we are building a container, this operation is fast, simply allocating a struct
1563  const RCP<panzer::ThyraObjContainer<Scalar> > thGlobalContainer =
1564  Teuchos::rcp_dynamic_cast<panzer::ThyraObjContainer<Scalar> >(ae_inargs.container_);
1565  const RCP<panzer::ThyraObjContainer<Scalar> > thGhostedContainer =
1566  Teuchos::rcp_dynamic_cast<panzer::ThyraObjContainer<Scalar> >(ae_inargs.ghostedContainer_);
1567 
1568  if (!Teuchos::is_null(f_out) && !Teuchos::is_null(W_out)) {
1569  PANZER_FUNC_TIME_MONITOR("panzer::ModelEvaluator::evalModel(f and J)");
1570 
1571  // only add auxiliary global data if Jacobian is being formed
1572  ae_inargs.addGlobalEvaluationData(nonParamGlobalEvaluationData_);
1573 
1574  // Set the targets
1575  thGlobalContainer->set_f_th(f_out);
1576  thGlobalContainer->set_A_th(W_out);
1577 
1578  // Zero values in ghosted container objects
1579  Thyra::assign(thGhostedContainer->get_f_th().ptr(),0.0);
1580  thGhostedContainer->initializeMatrix(0.0);
1581 
1582  ae_tm_.template getAsObject<panzer::Traits::Jacobian>()->evaluate(ae_inargs);
1583  }
1584  else if(!Teuchos::is_null(f_out) && Teuchos::is_null(W_out)) {
1585 
1586  PANZER_FUNC_TIME_MONITOR("panzer::ModelEvaluator::evalModel(f)");
1587 
1588  // don't add auxiliary global data if Jacobian is not computed.
1589  // this leads to zeroing of aux ops in special cases.
1590 
1591  thGlobalContainer->set_f_th(f_out);
1592 
1593  // Zero values in ghosted container objects
1594  Thyra::assign(thGhostedContainer->get_f_th().ptr(),0.0);
1595 
1596  ae_tm_.template getAsObject<panzer::Traits::Residual>()->evaluate(ae_inargs);
1597  }
1598  else if(Teuchos::is_null(f_out) && !Teuchos::is_null(W_out)) {
1599 
1600  PANZER_FUNC_TIME_MONITOR("panzer::ModelEvaluator::evalModel(J)");
1601 
1602  // only add auxiliary global data if Jacobian is being formed
1603  ae_inargs.addGlobalEvaluationData(nonParamGlobalEvaluationData_);
1604 
1605  // this dummy nonsense is needed only for scattering dirichlet conditions
1606  RCP<Thyra::VectorBase<Scalar> > dummy_f = Thyra::createMember(f_space_);
1607  thGlobalContainer->set_f_th(dummy_f);
1608  thGlobalContainer->set_A_th(W_out);
1609 
1610  // Zero values in ghosted container objects
1611  thGhostedContainer->initializeMatrix(0.0);
1612 
1613  ae_tm_.template getAsObject<panzer::Traits::Jacobian>()->evaluate(ae_inargs);
1614  }
1615 
1616  // HACK: set A to null before calling responses to avoid touching the
1617  // the Jacobian after it has been properly assembled. Should be fixed
1618  // by using a modified version of ae_inargs instead.
1619  thGlobalContainer->set_A_th(Teuchos::null);
1620 
1621  // Holding a rcp to f produces a seg fault in Rythmos when the next
1622  // f comes in and the resulting dtor is called. Need to discuss
1623  // with Ross. Clearing all references here works!
1624 
1625  thGlobalContainer->set_x_th(Teuchos::null);
1626  thGlobalContainer->set_dxdt_th(Teuchos::null);
1627  thGlobalContainer->set_f_th(Teuchos::null);
1628  thGlobalContainer->set_A_th(Teuchos::null);
1629 
1630  // reset parameters back to nominal values
1631  resetParameters();
1632 }
1633 
1634 template <typename Scalar>
1636 evalModelImpl_basic_g(const Thyra::ModelEvaluatorBase::InArgs<Scalar> &inArgs,
1637  const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs) const
1638 {
1639  // optional sanity check
1640  // TEUCHOS_ASSERT(required_basic_g(outArgs));
1641 
1642  // setup all the assembly in arguments (this is parameters and
1643  // x/x_dot). At this point with the exception of the one time dirichlet
1644  // beta that is all thats neccessary.
1645  panzer::AssemblyEngineInArgs ae_inargs;
1646  setupAssemblyInArgs(inArgs,ae_inargs);
1647 
1648  // set model parameters from supplied inArgs
1649  setParameters(inArgs);
1650 
1651  for(std::size_t i=0;i<responses_.size();i++) {
1652  Teuchos::RCP<Thyra::VectorBase<Scalar> > vec = outArgs.get_g(i);
1653  if(vec!=Teuchos::null) {
1654  std::string responseName = responses_[i]->name;
1656  = Teuchos::rcp_dynamic_cast<panzer::ResponseMESupportBase<panzer::Traits::Residual> >(
1657  responseLibrary_->getResponse<panzer::Traits::Residual>(responseName));
1658  resp->setVector(vec);
1659  }
1660  }
1661 
1662  // evaluator responses
1663  responseLibrary_->addResponsesToInArgs<panzer::Traits::Residual>(ae_inargs);
1664  responseLibrary_->evaluate<panzer::Traits::Residual>(ae_inargs);
1665 
1666  // reset parameters back to nominal values
1667  resetParameters();
1668 }
1669 
1670 template <typename Scalar>
1671 void
1673 evalModelImpl_basic_dgdx(const Thyra::ModelEvaluatorBase::InArgs<Scalar> &inArgs,
1674  const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs) const
1675 {
1676  typedef Thyra::ModelEvaluatorBase MEB;
1677 
1678  // optional sanity check
1679  TEUCHOS_ASSERT(required_basic_dgdx(outArgs));
1680 
1681  // set model parameters from supplied inArgs
1682  setParameters(inArgs);
1683 
1684  for(std::size_t i=0;i<responses_.size();i++) {
1685  // get "Vector" out of derivative, if its something else, throw an exception
1686  MEB::Derivative<Scalar> deriv = outArgs.get_DgDx(i);
1687  if(deriv.isEmpty())
1688  continue;
1689 
1690  Teuchos::RCP<Thyra::MultiVectorBase<Scalar> > vec = deriv.getMultiVector();
1691 
1692  if(vec!=Teuchos::null) {
1693 
1694  std::string responseName = responses_[i]->name;
1696  = Teuchos::rcp_dynamic_cast<panzer::ResponseMESupportBase<panzer::Traits::Jacobian> >(
1697  responseLibrary_->getResponse<panzer::Traits::Jacobian>(responseName));
1698  resp->setDerivative(vec);
1699  }
1700  }
1701 
1702  // setup all the assembly in arguments (this is parameters and
1703  // x/x_dot). At this point with the exception of the one time dirichlet
1704  // beta that is all thats neccessary.
1705  panzer::AssemblyEngineInArgs ae_inargs;
1706  setupAssemblyInArgs(inArgs,ae_inargs);
1707 
1708  // evaluate responses
1709  responseLibrary_->addResponsesToInArgs<panzer::Traits::Jacobian>(ae_inargs);
1710  responseLibrary_->evaluate<panzer::Traits::Jacobian>(ae_inargs);
1711 
1712  // reset parameters back to nominal values
1713  resetParameters();
1714 }
1715 
1716 template <typename Scalar>
1717 void
1719 evalModelImpl_basic_dgdp_scalar(const Thyra::ModelEvaluatorBase::InArgs<Scalar> &inArgs,
1720  const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs) const
1721 {
1722  using Teuchos::RCP;
1723  using Teuchos::rcp;
1724  using Teuchos::rcp_dynamic_cast;
1725 
1726  typedef Thyra::ModelEvaluatorBase MEB;
1727 
1728  // optional sanity check
1729  TEUCHOS_ASSERT(required_basic_dgdp_scalar(outArgs));
1730 
1731  // First find all of the active parameters for all responses
1732  std::vector<std::string> activeParameterNames;
1733  std::vector<int> activeParameters;
1734  int totalParameterCount = 0;
1735  for(std::size_t j=0; j<parameters_.size(); j++) {
1736 
1737  // skip non-scalar parameters
1738  if(parameters_[j]->is_distributed)
1739  continue;
1740 
1741  bool is_active = false;
1742  for(std::size_t i=0;i<responses_.size(); i++) {
1743 
1744  MEB::Derivative<Scalar> deriv = outArgs.get_DgDp(i,j);
1745  if(deriv.isEmpty())
1746  continue;
1747 
1748  Teuchos::RCP<Thyra::MultiVectorBase<Scalar> > vec = deriv.getMultiVector();
1749  if(vec!=Teuchos::null) {
1750  // get the response and tell it to fill the derivative vector
1751  std::string responseName = responses_[i]->name;
1754  responseLibrary_->getResponse<panzer::Traits::Tangent>(responseName));
1755  resp->setVector(vec);
1756  is_active = true;
1757  }
1758  }
1759 
1760  if (is_active) {
1761  for (std::size_t k=0; k<parameters_[j]->scalar_value.size(); k++) {
1762  std::string name = "PARAMETER_SENSITIVIES: "+(*parameters_[j]->names)[k];
1763  activeParameterNames.push_back(name);
1764  totalParameterCount++;
1765  }
1766  activeParameters.push_back(j);
1767  }
1768  }
1769 
1770  // setup all the assembly in arguments
1771  panzer::AssemblyEngineInArgs ae_inargs;
1772  setupAssemblyInArgs(inArgs,ae_inargs);
1773 
1774  // add active parameter names to assembly in-args
1775  RCP<panzer::GlobalEvaluationData> ged_activeParameters =
1776  rcp(new panzer::ParameterList_GlobalEvaluationData(activeParameterNames));
1777  ae_inargs.addGlobalEvaluationData("PARAMETER_NAMES",ged_activeParameters);
1778 
1779  // Initialize Fad components of all active parameters
1780  int paramIndex = 0;
1781  for (std::size_t ap=0; ap<activeParameters.size(); ++ap) {
1782  const int j = activeParameters[ap];
1783  for (unsigned int k=0; k < parameters_[j]->scalar_value.size(); k++) {
1784  panzer::Traits::FadType p(totalParameterCount, parameters_[j]->scalar_value[k].baseValue);
1785  p.fastAccessDx(paramIndex) = 1.0;
1786  parameters_[j]->scalar_value[k].family->template setValue<panzer::Traits::Tangent>(p);
1787  paramIndex++;
1788  }
1789  }
1790 
1791  // make sure that the total parameter count and the total parameter index match up
1792  TEUCHOS_ASSERT(paramIndex==totalParameterCount);
1793 
1794  // evaluate response tangent
1795  if(totalParameterCount>0) {
1796  responseLibrary_->addResponsesToInArgs<Traits::Tangent>(ae_inargs);
1797  responseLibrary_->evaluate<Traits::Tangent>(ae_inargs);
1798  }
1799 }
1800 
1801 template <typename Scalar>
1802 void
1804 evalModelImpl_basic_dgdp_distro(const Thyra::ModelEvaluatorBase::InArgs<Scalar> &inArgs,
1805  const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs) const
1806 {
1807  typedef Thyra::ModelEvaluatorBase MEB;
1808 
1809  // optional sanity check
1810  TEUCHOS_ASSERT(required_basic_dgdp_distro(outArgs));
1811 
1812  // loop over parameters, and then build a dfdp_rl only if they are distributed
1813  // and the user has provided the UGI. Note that this may be overly expensive if they
1814  // don't actually want those sensitivites because memory will be allocated unneccesarily.
1815  // It would be good to do this "just in time", but for now this is sufficient.
1816  for(std::size_t p=0;p<parameters_.size();p++) {
1817 
1818  // parameter is not distributed, a different path is
1819  // taken for those to compute dfdp
1820  if(!parameters_[p]->is_distributed)
1821  continue;
1822 
1823  ResponseLibrary<Traits> & rLibrary = *parameters_[p]->dgdp_rl;
1824 
1825  for(std::size_t r=0;r<responses_.size();r++) {
1826  // have derivatives been requested?
1827  MEB::Derivative<Scalar> deriv = outArgs.get_DgDp(r,p);
1828  if(deriv.isEmpty())
1829  continue;
1830 
1831  Teuchos::RCP<Thyra::MultiVectorBase<Scalar> > vec = deriv.getMultiVector();
1832 
1833  if(vec!=Teuchos::null) {
1834 
1835  // get the response and tell it to fill the derivative vector
1836  std::string responseName = responses_[r]->name;
1838  = Teuchos::rcp_dynamic_cast<panzer::ResponseMESupportBase<panzer::Traits::Jacobian> >(
1839  rLibrary.getResponse<panzer::Traits::Jacobian>(responseName));
1840 
1841  resp->setDerivative(vec);
1842  }
1843  }
1844 
1845  // setup all the assembly in arguments (this is parameters and x/x_dot).
1846  // make sure the correct seeding is performed
1847  panzer::AssemblyEngineInArgs ae_inargs;
1848  setupAssemblyInArgs(inArgs,ae_inargs);
1849 
1850  ae_inargs.first_sensitivities_name = (*parameters_[p]->names)[0]; // distributed parameters have one name!
1851  ae_inargs.gather_seeds.push_back(1.0); // this assumes that gather point is always the zero index of
1852  // gather seeds
1853 
1854  // evaluate responses
1855  rLibrary.addResponsesToInArgs<Traits::Jacobian>(ae_inargs);
1856  rLibrary.evaluate<Traits::Jacobian>(ae_inargs);
1857  }
1858 }
1859 
1860 template <typename Scalar>
1861 void
1863 evalModelImpl_basic_dfdp_scalar(const Thyra::ModelEvaluatorBase::InArgs<Scalar> &inArgs,
1864  const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs) const
1865 {
1866  using Teuchos::RCP;
1867  using Teuchos::rcp_dynamic_cast;
1868 
1869  typedef Thyra::ModelEvaluatorBase MEB;
1870 
1871  TEUCHOS_ASSERT(required_basic_dfdp_scalar(outArgs));
1872 
1873  // setup all the assembly in arguments (this is parameters and
1874  // x/x_dot). At this point with the exception of the one time dirichlet
1875  // beta that is all thats neccessary.
1876  panzer::AssemblyEngineInArgs ae_inargs;
1877  setupAssemblyInArgs(inArgs,ae_inargs);
1878 
1879  // First: Fill the output vectors from the input arguments structure. Put them
1880  // in the global evaluation data container so they are correctly communicated.
1882 
1883  std::vector<std::string> activeParameters;
1884 
1885  int totalParameterCount = 0;
1886  for(std::size_t i=0; i < parameters_.size(); i++) {
1887  // skip non-scalar parameters
1888  if(parameters_[i]->is_distributed)
1889  continue;
1890 
1891  // have derivatives been requested?
1892  MEB::Derivative<Scalar> deriv = outArgs.get_DfDp(i);
1893  if(deriv.isEmpty())
1894  continue;
1895 
1896  // grab multivector, make sure its the right dimension
1897  Teuchos::RCP<Thyra::MultiVectorBase<Scalar> > mVec = deriv.getMultiVector();
1898  TEUCHOS_ASSERT(mVec->domain()->dim()==Teuchos::as<int>(parameters_[i]->scalar_value.size()));
1899 
1900  for (std::size_t j=0; j < parameters_[i]->scalar_value.size(); j++) {
1901 
1902  // build containers for each vector
1905  RCP<LinearObjContainer> globalContainer = loc_pair->getGlobalLOC();
1906 
1907  // stuff target vector into global container
1908  RCP<Thyra::VectorBase<Scalar> > vec = mVec->col(j);
1909  RCP<panzer::ThyraObjContainer<Scalar> > thGlobalContainer =
1910  Teuchos::rcp_dynamic_cast<panzer::ThyraObjContainer<Scalar> >(globalContainer);
1911  thGlobalContainer->set_f_th(vec);
1912 
1913  // add container into in args object
1914  std::string name = "PARAMETER_SENSITIVIES: "+(*parameters_[i]->names)[j];
1915  ae_inargs.addGlobalEvaluationData(name,loc_pair->getGhostedLOC());
1916  ae_inargs.addGlobalEvaluationData(name+"_pair",loc_pair);
1917 
1918  activeParameters.push_back(name);
1919  totalParameterCount++;
1920  }
1921  }
1922 
1923  // Second: For all parameters that require derivative sensitivities, put in a name
1924  // so that the scatter can realize which sensitivity vectors it needs to fill
1926 
1927  RCP<GlobalEvaluationData> ged_activeParameters
1928  = Teuchos::rcp(new ParameterList_GlobalEvaluationData(activeParameters));
1929  ae_inargs.addGlobalEvaluationData("PARAMETER_NAMES",ged_activeParameters);
1930 
1931  // Third: Now seed all the parameters in the parameter vector so that derivatives
1932  // can be properly computed.
1934 
1935  int paramIndex = 0;
1936  for(std::size_t i=0; i < parameters_.size(); i++) {
1937  // skip non-scalar parameters
1938  if(parameters_[i]->is_distributed)
1939  continue;
1940 
1941  // don't modify the parameter if its not needed
1942  MEB::Derivative<Scalar> deriv = outArgs.get_DfDp(i);
1943  if(deriv.isEmpty()) {
1944  // reinitialize values that should not have sensitivities computed (this is a precaution)
1945  for (unsigned int j=0; j < parameters_[i]->scalar_value.size(); j++) {
1946  Traits::FadType p = Traits::FadType(totalParameterCount,
1947  parameters_[i]->scalar_value[j].baseValue);
1948  parameters_[i]->scalar_value[j].family->template setValue<panzer::Traits::Tangent>(p);
1949  }
1950  continue;
1951  }
1952  else {
1953  // loop over each parameter in the vector, initializing the AD type
1954  for (unsigned int j=0; j < parameters_[i]->scalar_value.size(); j++) {
1955  Traits::FadType p = Traits::FadType(totalParameterCount,
1956  parameters_[i]->scalar_value[j].baseValue);
1957  p.fastAccessDx(paramIndex) = 1.0;
1958  parameters_[i]->scalar_value[j].family->template setValue<panzer::Traits::Tangent>(p);
1959  paramIndex++;
1960  }
1961  }
1962  }
1963 
1964  // make sure that the total parameter count and the total parameter index match up
1965  TEUCHOS_ASSERT(paramIndex==totalParameterCount);
1966 
1967  // Fourth: Actually evaluate the residual's sensitivity to the parameters
1969 
1970  if(totalParameterCount>0) {
1971  PANZER_FUNC_TIME_MONITOR("panzer::ModelEvaluator::evalModel(df/dp)");
1972  ae_tm_.getAsObject<panzer::Traits::Tangent>()->evaluate(ae_inargs);
1973  }
1974 }
1975 
1976 template <typename Scalar>
1977 void
1979 evalModelImpl_basic_dfdp_scalar_fd(const Thyra::ModelEvaluatorBase::InArgs<Scalar> &inArgs,
1980  const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs) const
1981 {
1982  PANZER_FUNC_TIME_MONITOR("panzer::ModelEvaluator::evalModel(df/dp)");
1983 
1984  using Teuchos::RCP;
1985  using Teuchos::rcp_dynamic_cast;
1986 
1987  typedef Thyra::ModelEvaluatorBase MEB;
1988 
1989  TEUCHOS_ASSERT(required_basic_dfdp_scalar(outArgs));
1990 
1991  // First evaluate the model without df/dp for the base point
1992  // Maybe it would be better to set all outArgs and then remove the df/dp ones,
1993  // but I couldn't get that to work.
1994  MEB::OutArgs<Scalar> outArgs_base = this->createOutArgs();
1995  if (outArgs.get_f() == Teuchos::null)
1996  outArgs_base.set_f(Thyra::createMember(this->get_f_space()));
1997  else
1998  outArgs_base.set_f(outArgs.get_f());
1999  outArgs_base.set_W_op(outArgs.get_W_op());
2000  this->evalModel(inArgs, outArgs_base);
2001  RCP<const Thyra::VectorBase<Scalar> > f = outArgs_base.get_f();
2002  RCP<const Thyra::VectorBase<Scalar> > x = inArgs.get_x();
2004  if (inArgs.supports(MEB::IN_ARG_x_dot))
2005  x_dot = inArgs.get_x_dot();
2006 
2007  // Create in/out args for FD calculation
2008  RCP<Thyra::VectorBase<Scalar> > fd = Thyra::createMember(this->get_f_space());
2009  MEB::OutArgs<Scalar> outArgs_fd = this->createOutArgs();
2010  outArgs_fd.set_f(fd);
2011 
2012  RCP<Thyra::VectorBase<Scalar> > xd = Thyra::createMember(this->get_x_space());
2014  if (x_dot != Teuchos::null)
2015  xd_dot = Thyra::createMember(this->get_x_space());
2016  MEB::InArgs<Scalar> inArgs_fd = this->createInArgs();
2017  inArgs_fd.setArgs(inArgs); // This sets all inArgs that we don't override below
2018  inArgs_fd.set_x(xd);
2019  if (x_dot != Teuchos::null)
2020  inArgs_fd.set_x_dot(xd_dot);
2021 
2022  const double h = fd_perturb_size_;
2023  for(std::size_t i=0; i < parameters_.size(); i++) {
2024 
2025  // skip non-scalar parameters
2026  if(parameters_[i]->is_distributed)
2027  continue;
2028 
2029  // have derivatives been requested?
2030  MEB::Derivative<Scalar> deriv = outArgs.get_DfDp(i);
2031  if(deriv.isEmpty())
2032  continue;
2033 
2034  // grab multivector, make sure its the right dimension
2035  RCP<Thyra::MultiVectorBase<Scalar> > dfdp = deriv.getMultiVector();
2036  TEUCHOS_ASSERT(dfdp->domain()->dim()==Teuchos::as<int>(parameters_[i]->scalar_value.size()));
2037 
2038  // Get parameter vector and tangent vectors
2039  RCP<const Thyra::VectorBase<Scalar> > p = inArgs.get_p(i);
2040  RCP<const Thyra::VectorBase<Scalar> > dx_v = inArgs.get_p(i+parameters_.size());
2042  rcp_dynamic_cast<const Thyra::DefaultMultiVectorProductVector<Scalar> >(dx_v,true)->getMultiVector();
2045  if (x_dot != Teuchos::null) {
2046  dx_dot_v =inArgs.get_p(i+parameters_.size()+tangent_space_.size());
2047  dx_dot =
2048  rcp_dynamic_cast<const Thyra::DefaultMultiVectorProductVector<Scalar> >(dx_dot_v,true)->getMultiVector();
2049  }
2050 
2051  // Create perturbed parameter vector
2052  RCP<Thyra::VectorBase<Scalar> > pd = Thyra::createMember(this->get_p_space(i));
2053  inArgs_fd.set_p(i,pd);
2054 
2055  for (std::size_t j=0; j < parameters_[i]->scalar_value.size(); j++) {
2056 
2057  // Perturb parameter vector
2058  Thyra::copy(*p, pd.ptr());
2059  Thyra::set_ele(j, Thyra::get_ele(*p,j)+h, pd.ptr());
2060 
2061  // Perturb state vectors using tangents
2062  Thyra::V_VpStV(xd.ptr(), *x, h, *(dx)->col(j));
2063  if (x_dot != Teuchos::null)
2064  Thyra::V_VpStV(xd_dot.ptr(), *x_dot, h, *(dx_dot)->col(j));
2065 
2066  // Evaluate perturbed residual
2067  Thyra::assign(fd.ptr(), 0.0);
2068  this->evalModel(inArgs_fd, outArgs_fd);
2069 
2070  // FD calculation
2071  Thyra::V_StVpStV(dfdp->col(j).ptr(), 1.0/h, *fd, -1.0/h, *f);
2072 
2073  // Reset parameter back to un-perturbed value
2074  parameters_[i]->scalar_value[j].family->setRealValueForAllTypes(Thyra::get_ele(*p,j));
2075 
2076  }
2077  }
2078 }
2079 
2080 template <typename Scalar>
2081 void
2083 evalModelImpl_basic_dfdp_distro(const Thyra::ModelEvaluatorBase::InArgs<Scalar> &inArgs,
2084  const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs) const
2085 {
2086  using Teuchos::RCP;
2087  using Teuchos::rcp_dynamic_cast;
2088  using Teuchos::null;
2089 
2090  typedef Thyra::ModelEvaluatorBase MEB;
2091 
2092  TEUCHOS_ASSERT(required_basic_dfdp_distro(outArgs));
2093 
2094  // loop over parameters, and then build a dfdp_rl only if they are distributed
2095  // and the user has provided the UGI. Note that this may be overly expensive if they
2096  // don't actually want those sensitivites because memory will be allocated unneccesarily.
2097  // It would be good to do this "just in time", but for now this is sufficient.
2098  for(std::size_t p=0;p<parameters_.size();p++) {
2099 
2100  // parameter is not distributed, a different path is
2101  // taken for those to compute dfdp
2102  if(!parameters_[p]->is_distributed)
2103  continue;
2104 
2105  // parameter is distributed but has no global indexer.
2106  // thus the user doesn't want sensitivities!
2107  if(parameters_[p]->dfdp_rl==null)
2108  continue;
2109 
2110  // have derivatives been requested?
2111  MEB::Derivative<Scalar> deriv = outArgs.get_DfDp(p);
2112  if(deriv.isEmpty())
2113  continue;
2114 
2115  ResponseLibrary<Traits> & rLibrary = *parameters_[p]->dfdp_rl;
2116 
2117  // get the response and tell it to fill the derivative operator
2118  RCP<Response_Residual<Traits::Jacobian> > response_jacobian =
2119  rcp_dynamic_cast<Response_Residual<Traits::Jacobian> >(rLibrary.getResponse<Traits::Jacobian>("RESIDUAL"));
2120  response_jacobian->setJacobian(deriv.getLinearOp());
2121 
2122  // setup all the assembly in arguments (this is parameters and x/x_dot).
2123  // make sure the correct seeding is performed
2124  panzer::AssemblyEngineInArgs ae_inargs;
2125  setupAssemblyInArgs(inArgs,ae_inargs);
2126 
2127  ae_inargs.first_sensitivities_name = (*parameters_[p]->names)[0]; // distributed parameters have one name!
2128  ae_inargs.gather_seeds.push_back(1.0); // this assumes that gather point is always the zero index of
2129  // gather seeds
2130  rLibrary.addResponsesToInArgs<Traits::Jacobian>(ae_inargs);
2131 
2132  rLibrary.evaluate<Traits::Jacobian>(ae_inargs);
2133  }
2134 }
2135 
2136 template <typename Scalar>
2138 required_basic_g(const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs) const
2139 {
2140  // determine if any of the outArgs are not null!
2141  bool activeGArgs = false;
2142  for(int i=0;i<outArgs.Ng();i++)
2143  activeGArgs |= (outArgs.get_g(i)!=Teuchos::null);
2144 
2145  return activeGArgs | required_basic_dgdx(outArgs);
2146 }
2147 
2148 template <typename Scalar>
2150 required_basic_dgdx(const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs) const
2151 {
2152  typedef Thyra::ModelEvaluatorBase MEB;
2153 
2154  // determine if any of the outArgs are not null!
2155  bool activeGArgs = false;
2156  for(int i=0;i<outArgs.Ng();i++) {
2157  // no derivatives are supported
2158  if(outArgs.supports(MEB::OUT_ARG_DgDx,i).none())
2159  continue;
2160 
2161  // this is basically a redundant computation
2162  activeGArgs |= (!outArgs.get_DgDx(i).isEmpty());
2163  }
2164 
2165  return activeGArgs;
2166 }
2167 
2168 template <typename Scalar>
2170 required_basic_dgdp_scalar(const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs) const
2171 {
2172  typedef Thyra::ModelEvaluatorBase MEB;
2173 
2174  // determine if any of the outArgs are not null!
2175  bool activeGArgs = false;
2176  for(int i=0;i<outArgs.Ng();i++) {
2177  for(int p=0;p<Teuchos::as<int>(parameters_.size());p++) {
2178 
2179  // only look at scalar parameters
2180  if(parameters_[p]->is_distributed)
2181  continue;
2182 
2183  // no derivatives are supported
2184  if(outArgs.supports(MEB::OUT_ARG_DgDp,i,p).none())
2185  continue;
2186 
2187  activeGArgs |= (!outArgs.get_DgDp(i,p).isEmpty());
2188  }
2189  }
2190 
2191  return activeGArgs;
2192 }
2193 
2194 template <typename Scalar>
2196 required_basic_dgdp_distro(const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs) const
2197 {
2198  typedef Thyra::ModelEvaluatorBase MEB;
2199 
2200  // determine if any of the outArgs are not null!
2201  bool activeGArgs = false;
2202  for(int i=0;i<outArgs.Ng();i++) {
2203  for(int p=0;p<Teuchos::as<int>(parameters_.size());p++) {
2204 
2205  // only look at distributed parameters
2206  if(!parameters_[p]->is_distributed)
2207  continue;
2208 
2209  // no derivatives are supported
2210  if(outArgs.supports(MEB::OUT_ARG_DgDp,i,p).none())
2211  continue;
2212 
2213  activeGArgs |= (!outArgs.get_DgDp(i,p).isEmpty());
2214  }
2215  }
2216 
2217  return activeGArgs;
2218 }
2219 
2220 template <typename Scalar>
2222 required_basic_dfdp_scalar(const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs) const
2223 {
2224  typedef Thyra::ModelEvaluatorBase MEB;
2225 
2226  // determine if any of the outArgs are not null!
2227  bool activeFPArgs = false;
2228  for(int i=0;i<Teuchos::as<int>(parameters_.size());i++) {
2229 
2230  // this is for scalar parameters only
2231  if(parameters_[i]->is_distributed)
2232  continue;
2233 
2234  // no derivatives are supported
2235  if(outArgs.supports(MEB::OUT_ARG_DfDp,i).none())
2236  continue;
2237 
2238  // this is basically a redundant computation
2239  activeFPArgs |= (!outArgs.get_DfDp(i).isEmpty());
2240  }
2241 
2242  return activeFPArgs;
2243 }
2244 
2245 template <typename Scalar>
2247 required_basic_dfdp_distro(const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs) const
2248 {
2249  typedef Thyra::ModelEvaluatorBase MEB;
2250 
2251  // determine if any of the outArgs are not null!
2252  bool activeFPArgs = false;
2253  for(int i=0;i<Teuchos::as<int>(parameters_.size());i++) {
2254 
2255  // this is for scalar parameters only
2256  if(!parameters_[i]->is_distributed)
2257  continue;
2258 
2259  // no derivatives are supported
2260  if(outArgs.supports(MEB::OUT_ARG_DfDp,i).none())
2261  continue;
2262 
2263  // this is basically a redundant computation
2264  activeFPArgs |= (!outArgs.get_DfDp(i).isEmpty());
2265  }
2266 
2267  return activeFPArgs;
2268 }
2269 
2270 template <typename Scalar>
2274  const std::vector<Teuchos::RCP<panzer::PhysicsBlock> >& physicsBlocks,
2275  const std::vector<panzer::BC> & bcs,
2276  const panzer::EquationSetFactory & eqset_factory,
2277  const panzer::BCStrategyFactory& bc_factory,
2279  const Teuchos::ParameterList& closure_models,
2280  const Teuchos::ParameterList& user_data,
2281  const bool write_graphviz_file,
2282  const std::string& graphviz_file_prefix)
2283 {
2284  using Teuchos::RCP;
2285  using Teuchos::rcp;
2286  using Teuchos::null;
2287 
2288  // loop over parameters, and then build a dfdp_rl only if they are distributed
2289  // and the user has provided the UGI. Note that this may be overly expensive if they
2290  // don't actually want those sensitivites because memory will be allocated unneccesarily.
2291  // It would be good to do this "just in time", but for now this is sufficient.
2292  for(std::size_t p=0;p<parameters_.size();p++) {
2293  // parameter is not distributed, a different path is
2294  // taken for those to compute dfdp
2295  if(!parameters_[p]->is_distributed)
2296  continue;
2297 
2298  // parameter is distributed but has no global indexer.
2299  // thus the user doesn't want sensitivities!
2300  if(parameters_[p]->global_indexer==null)
2301  continue;
2302 
2303  // build the linear object factory that has the correct sizing for
2304  // the sensitivity matrix (parameter sized domain, residual sized range)
2306  parameters_[p]->global_indexer);
2307 
2308  // the user wants global sensitivities, hooray! Build and setup the response library
2309  RCP<ResponseLibrary<Traits> > rLibrary
2310  = Teuchos::rcp(new ResponseLibrary<Traits>(wc,lof_->getRangeGlobalIndexer(),
2311  param_lof,true));
2312  rLibrary->buildResidualResponseEvaluators(physicsBlocks,eqset_factory,bcs,bc_factory,
2313  cm_factory,closure_models,user_data,
2314  write_graphviz_file,graphviz_file_prefix);
2315 
2316  // make sure parameter response library is correct
2317  parameters_[p]->dfdp_rl = rLibrary;
2318  }
2319 }
2320 
2321 template <typename Scalar>
2325  const std::vector<Teuchos::RCP<panzer::PhysicsBlock> >& physicsBlocks,
2326  const std::vector<panzer::BC>& /* bcs */,
2327  const panzer::EquationSetFactory & eqset_factory,
2328  const panzer::BCStrategyFactory& /* bc_factory */,
2330  const Teuchos::ParameterList& closure_models,
2331  const Teuchos::ParameterList& user_data,
2332  const bool write_graphviz_file,
2333  const std::string& graphviz_file_prefix)
2334 {
2335  using Teuchos::RCP;
2336  using Teuchos::rcp;
2337  using Teuchos::null;
2338 
2339  // loop over parameters, and then build a dfdp_rl only if they are distributed
2340  // and the user has provided the UGI. Note that this may be overly expensive if they
2341  // don't actually want those sensitivites because memory will be allocated unneccesarily.
2342  // It would be good to do this "just in time", but for now this is sufficient.
2343  for(std::size_t p=0;p<parameters_.size();p++) {
2344  // parameter is not distributed, a different path is
2345  // taken for those to compute dfdp
2346  if(!parameters_[p]->is_distributed)
2347  continue;
2348 
2349  // parameter is distributed but has no global indexer.
2350  // thus the user doesn't want sensitivities!
2351  if(parameters_[p]->global_indexer==null)
2352  continue;
2353 
2354  // extract the linear object factory that has the correct sizing for
2355  // the sensitivity vector
2356  RCP<const LinearObjFactory<Traits> > param_lof = parameters_[p]->dfdp_rl->getLinearObjFactory();
2357  RCP<const GlobalIndexer > param_ugi = parameters_[p]->global_indexer;
2358 
2359  // the user wants global sensitivities, hooray! Build and setup the response library
2360  RCP<ResponseLibrary<Traits> > rLibrary
2361  = Teuchos::rcp(new ResponseLibrary<Traits>(wc,lof_->getRangeGlobalIndexer(), lof_));
2362 
2363 
2364  // build evaluators for all flexible responses
2365  for(std::size_t r=0;r<responses_.size();r++) {
2366  // only responses with a builder are non null!
2367  if(responses_[r]->builder==Teuchos::null)
2368  continue;
2369 
2370  // set the current derivative information in the builder
2371  // responses_[r]->builder->setDerivativeInformationBase(param_lof,param_ugi);
2372  responses_[r]->builder->setDerivativeInformation(param_lof);
2373 
2374  // add the response
2375  rLibrary->addResponse(responses_[r]->name,
2376  responses_[r]->wkst_desc,
2377  *responses_[r]->builder);
2378  }
2379 
2380  rLibrary->buildResponseEvaluators(physicsBlocks,eqset_factory,
2381  cm_factory,closure_models,user_data,
2382  write_graphviz_file,graphviz_file_prefix);
2383 
2384  // make sure parameter response library is correct
2385  parameters_[p]->dgdp_rl = rLibrary;
2386  }
2387 }
2388 
2389 template <typename Scalar>
2391 setOneTimeDirichletBeta(const Scalar & beta) const
2392 {
2393  oneTimeDirichletBeta_on_ = true;
2394  oneTimeDirichletBeta_ = beta;
2395 }
2396 
2397 template <typename Scalar>
2401  const Teuchos::Array<Scalar> & in_values) const
2402 {
2403  using Teuchos::RCP;
2404  using Teuchos::rcp;
2405  using Teuchos::rcp_dynamic_cast;
2406  using Teuchos::ptrFromRef;
2407 
2408  TEUCHOS_ASSERT(in_names.size()==in_values.size());
2409 
2410  // Check that the parameters are valid (i.e., they already exist in the parameter library)
2411  // std::size_t np = in_names.size();
2412  // for(std::size_t i=0;i<np;i++)
2413  // TEUCHOS_TEST_FOR_EXCEPTION(!global_data_->pl->isParameter(in_names[i]),
2414  // std::logic_error,
2415  // "Parameter \"" << in_names[i] << "\" does not exist in parameter library!");
2416 
2417  RCP<ParameterObject> paramObj = rcp(new ParameterObject);
2418 
2419  paramObj->names = rcp(new Teuchos::Array<std::string>(in_names));
2420  paramObj->is_distributed = false;
2421 
2422  // register all the scalar parameters, setting initial
2423  for(int i=0;i<in_names.size();i++)
2424  registerScalarParameter(in_names[i],*global_data_->pl,in_values[i]);
2425 
2426  paramObj->scalar_value = panzer::ParamVec();
2427  global_data_->pl->fillVector<panzer::Traits::Residual>(*paramObj->names, paramObj->scalar_value);
2428 
2429  // build initial condition vector
2430  paramObj->space =
2431  Thyra::locallyReplicatedDefaultSpmdVectorSpace<Scalar>(
2432  rcp(new Teuchos::MpiComm<long int>(lof_->getComm().getRawMpiComm())),paramObj->names->size());
2433 
2434  // fill vector with parameter values
2436  RCP<Thyra::VectorBase<Scalar> > initial_value = Thyra::createMember(paramObj->space);
2437  RCP<Thyra::SpmdVectorBase<Scalar> > vec = rcp_dynamic_cast<Thyra::SpmdVectorBase<Scalar> >(initial_value);
2438  vec->getNonconstLocalData(ptrFromRef(data));
2439  for (unsigned int i=0; i < paramObj->scalar_value.size(); i++)
2440  data[i] = in_values[i];
2441 
2442  paramObj->initial_value = initial_value;
2443 
2444  return paramObj;
2445 }
2446 
2447 template <typename Scalar>
2450 createDistributedParameter(const std::string & key,
2451  const Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> > & vs,
2452  const Teuchos::RCP<const Thyra::VectorBase<Scalar> > & initial,
2453  const Teuchos::RCP<const GlobalIndexer> & ugi) const
2454 {
2455  using Teuchos::RCP;
2456  using Teuchos::rcp;
2457 
2458  RCP<ParameterObject> paramObj = rcp(new ParameterObject);
2459 
2460  paramObj->is_distributed = true;
2461  paramObj->names = rcp(new Teuchos::Array<std::string>());
2462  paramObj->names->push_back(key);
2463  paramObj->space = vs;
2464  paramObj->initial_value = initial;
2465 
2466  paramObj->global_indexer = ugi;
2467 
2468  return paramObj;
2469 }
2470 
2471 template <typename Scalar>
2472 void
2474 setParameters(const Thyra::ModelEvaluatorBase::InArgs<Scalar> &inArgs) const
2475 {
2476  for(std::size_t i=0; i < parameters_.size(); i++) {
2477 
2478  // skip non-scalar parameters (for now)
2479  if(parameters_[i]->is_distributed)
2480  continue;
2481 
2482  // set parameter values for given parameter vector for all evaluation types
2483  Teuchos::RCP<const Thyra::VectorBase<Scalar> > p = inArgs.get_p(i);
2484  if (p != Teuchos::null) {
2485  for (unsigned int j=0; j < parameters_[i]->scalar_value.size(); j++) {
2486  parameters_[i]->scalar_value[j].family->setRealValueForAllTypes(Thyra::get_ele(*p,j));
2487  }
2488  }
2489 
2490  }
2491 }
2492 
2493 template <typename Scalar>
2494 void
2497 {
2498  for(std::size_t i=0; i < parameters_.size(); i++) {
2499 
2500  // skip non-scalar parameters (for now)
2501  if(parameters_[i]->is_distributed)
2502  continue;
2503 
2504  // Reset each parameter back to its nominal
2505  for (unsigned int j=0; j < parameters_[i]->scalar_value.size(); j++) {
2506  parameters_[i]->scalar_value[j].family->setRealValueForAllTypes(Thyra::get_ele(*(parameters_[i]->initial_value),j));
2507  }
2508 
2509  }
2510 }
2511 
2512 #endif // __Panzer_ModelEvaluator_impl_hpp__
Teuchos::RCP< const LinearObjFactory< panzer::Traits > > cloneWithNewDomain(const LinearObjFactory< panzer::Traits > &lof, const Teuchos::RCP< const GlobalIndexer > &dUgi)
Clone a linear object factory, but using a different domain.
Interface for constructing a BCStrategy_TemplateManager.
void setupVolumeFieldManagers(const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &physicsBlocks, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &cm_factory, const Teuchos::ParameterList &closure_models, const LinearObjFactory< panzer::Traits > &lo_factory, const Teuchos::ParameterList &user_data)
Allocates and initializes an equation set template manager.
virtual void evalModelImpl_basic_dgdp_distro(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int i) const override
Teuchos::RCP< panzer::ResponseLibrary< panzer::Traits > > dfdp_rl
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_x_space() const override
Teuchos::RCP< const GlobalIndexer > global_indexer
void evalModel_D2gDp2(int rIndex, int pIndex, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &delta_x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &D2gDp2) const
void initializeNominalValues() const
Initialize the nominal values with good starting conditions.
int addDistributedParameter(const std::string &name, const Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > &vs, const Teuchos::RCP< GlobalEvaluationData > &ged, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &initial, const Teuchos::RCP< const GlobalIndexer > &ugi=Teuchos::null)
void addResponsesToInArgs(panzer::AssemblyEngineInArgs &input_args) const
int addFlexibleResponse(const std::string &responseName, const std::vector< WorksetDescriptor > &wkst_desc, const Teuchos::RCP< ResponseMESupportBuilderBase > &builder)
T & get(const std::string &name, T def_value)
void evalModel_D2gDx2(int rIndex, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &delta_x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &D2gDx2) const
bool is_null(const std::shared_ptr< T > &p)
void evalModel_D2fDp2(int pIndex, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &delta_x, const Teuchos::RCP< Thyra::LinearOpBase< Scalar > > &D2fDp2) const
void writeVolumeGraphvizDependencyFiles(std::string filename_prefix, const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &physicsBlocks) const
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_f_space() const override
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< ParameterObject > createScalarParameter(const Teuchos::Array< std::string > &names, const Teuchos::Array< Scalar > &in_values) const
panzer::AssemblyEngine_TemplateManager< panzer::Traits > ae_tm_
void buildDistroParamDgDp_RL(const Teuchos::RCP< panzer::WorksetContainer > &wc, const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &physicsBlocks, const std::vector< panzer::BC > &bcs, const panzer::EquationSetFactory &eqset_factory, const panzer::BCStrategyFactory &bc_factory, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &cm_factory, const Teuchos::ParameterList &closure_models, const Teuchos::ParameterList &user_data, const bool write_graphviz_file=false, const std::string &graphviz_file_prefix="")
Teuchos::RCP< LinearObjContainer > getGlobalLOC() const
Teuchos::RCP< Thyra::LinearOpBase< Scalar > > create_DfDp_op(int i) const override
void evalModel_D2fDxDp(int pIndex, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &delta_p, const Teuchos::RCP< Thyra::LinearOpBase< Scalar > > &D2fDxDp) const
void setupAssemblyInArgs(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, panzer::AssemblyEngineInArgs &ae_inargs) const
void buildVolumeFieldManagers(const bool value)
Thyra::ModelEvaluatorBase::InArgs< Scalar > createInArgs() const override
void writeBCGraphvizDependencyFiles(std::string filename_prefix) const
void evalModel_D2gDpDx(int rIndex, int pIndex, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &delta_x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &D2gDpDx) const
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_g_space(int i) const override
void setOneTimeDirichletBeta(const Scalar &beta) const
virtual void evalModelImpl_basic(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
Evaluate a simple model, meaning a residual and a jacobian, no fancy stochastic galerkin or multipoin...
void evalModel_D2fDx2(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &delta_x, const Teuchos::RCP< Thyra::LinearOpBase< Scalar > > &D2fDx2) const
bool required_basic_g(const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
Does this set of out args require a simple response?
Sacado::ScalarParameterVector< panzer::EvaluationTraits > ParamVec
Teuchos::RCP< panzer::LinearObjContainer > ghostedContainer_
virtual void evalModelImpl_basic_dgdx(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
bool required_basic_dgdp_distro(const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
Are their required responses in the out args? DgDp.
Teuchos::RCP< Epetra_MultiVector > getMultiVector(const std::string &modelEvalDescription, const ModelEvaluator::Derivative &deriv, const std::string &derivName, const ModelEvaluator::EDerivativeMultiVectorOrientation mvOrientation)
const std::string & get_g_name(int i) const
Teuchos::RCP< Thyra::LinearOpBase< Scalar > > create_W_op() const override
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > x_space_
void setParameters(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs) const
void evaluate(const panzer::AssemblyEngineInArgs &input_args)
Teuchos::RCP< LinearObjContainer > getGhostedLOC() const
Teuchos::RCP< panzer::LinearObjContainer > container_
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > f_space_
Teuchos::ArrayView< const std::string > get_g_names(int i) const override
Thyra::ModelEvaluatorBase::OutArgs< Scalar > createOutArgsImpl() const override
void addGlobalEvaluationData(const std::string &key, const Teuchos::RCP< GlobalEvaluationData > &ged)
void setWorksetContainer(const Teuchos::RCP< WorksetContainer > &wc)
bool isParameter(const std::string &name) const
Teuchos::RCP< const Thyra::LinearOpWithSolveFactoryBase< Scalar > > get_W_factory() const override
PHX::MDField< ScalarT, panzer::Cell, panzer::IP > result
A field that will be used to build up the result of the integral we&#39;re performing.
virtual void evalModelImpl_basic_g(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
Construct a simple response dicatated by this set of out args.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< ParameterObject > createDistributedParameter(const std::string &key, const Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > &vs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &initial, const Teuchos::RCP< const GlobalIndexer > &ugi) const
Teuchos::RCP< const panzer::LinearObjFactory< panzer::Traits > > lof_
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_p_space(int i) const override
bool required_basic_dgdp_scalar(const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
Are their required responses in the out args? DgDp.
Ptr< T > ptr() const
bool required_basic_dgdx(const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
Are their required responses in the out args? DgDx.
virtual void evalModelImpl(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const override
virtual void setDerivativeInformation(const Teuchos::RCP< const panzer::LinearObjFactory< panzer::Traits > > &linearObjFactory)=0
Teuchos::RCP< ResponseBase > getResponse(const std::string &responseName) const
virtual void evalModelImpl_basic_dfdp_distro(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
virtual void evalModelImpl_basic_dfdp_scalar(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
void buildDistroParamDfDp_RL(const Teuchos::RCP< panzer::WorksetContainer > &wc, const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &physicsBlocks, const std::vector< panzer::BC > &bcs, const panzer::EquationSetFactory &eqset_factory, const panzer::BCStrategyFactory &bc_factory, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &cm_factory, const Teuchos::ParameterList &closure_models, const Teuchos::ParameterList &user_data, const bool write_graphviz_file=false, const std::string &graphviz_file_prefix="")
void setupModel(const Teuchos::RCP< panzer::WorksetContainer > &wc, const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &physicsBlocks, const std::vector< panzer::BC > &bcs, const panzer::EquationSetFactory &eqset_factory, const panzer::BCStrategyFactory &bc_factory, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &volume_cm_factory, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &bc_cm_factory, const Teuchos::ParameterList &closure_models, const Teuchos::ParameterList &user_data, bool writeGraph=false, const std::string &graphPrefix="", const Teuchos::ParameterList &me_params=Teuchos::ParameterList())
virtual void evalModelImpl_basic_dfdp_scalar_fd(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
PANZER_FADTYPE FadType
Thyra::ModelEvaluatorBase::InArgs< Scalar > getNominalValues() const override
void push_back(const value_type &x)
size_type size() const
Teuchos::RCP< panzer::ResponseLibrary< panzer::Traits > > responseLibrary_
void addNonParameterGlobalEvaluationData(const std::string &name, const Teuchos::RCP< GlobalEvaluationData > &ged)
bool required_basic_dfdp_scalar(const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
Are derivatives of the residual with respect to the scalar parameters in the out args? DfDp.
void evalModel_D2fDpDx(int pIndex, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &delta_x, const Teuchos::RCP< Thyra::LinearOpBase< Scalar > > &D2fDpDx) const
void evalModel_D2gDxDp(int rIndex, int pIndex, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &delta_p, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &D2gDxDp) const
void setupBCFieldManagers(const std::vector< panzer::BC > &bcs, const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &physicsBlocks, const panzer::EquationSetFactory &eqset_factory, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &cm_factory, const panzer::BCStrategyFactory &bc_factory, const Teuchos::ParameterList &closure_models, const LinearObjFactory< panzer::Traits > &lo_factory, const Teuchos::ParameterList &user_data)
void applyDirichletBCs(const Teuchos::RCP< Thyra::VectorBase< Scalar > > &x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &f) const
#define TEUCHOS_ASSERT(assertion_test)
bool required_basic_dfdp_distro(const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
Are derivatives of the residual with respect to the distributed parameters in the out args...
void registerScalarParameter(const std::string name, panzer::ParamLib &pl, double realValue)
void setActiveEvaluationTypes(const std::vector< bool > &aet)
Set a vector of active evaluation types to allocate.
virtual void evalModelImpl_basic_dgdp_scalar(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
void buildBCFieldManagers(const bool value)
int addParameter(const std::string &name, const Scalar &initial)