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