Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
CDR_Model_impl.hpp
Go to the documentation of this file.
1 //@HEADER
2 // *****************************************************************************
3 // Tempus: Time Integration and Sensitivity Analysis Package
4 //
5 // Copyright 2017 NTESS and the Tempus contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 //@HEADER
9 
10 #ifndef NOX_THYRA_MODEL_EVALUATOR_1DFEM_DEF_HPP
11 #define NOX_THYRA_MODEL_EVALUATOR_1DFEM_DEF_HPP
12 
13 // Thyra support
14 #include "Thyra_DefaultSpmdVectorSpace.hpp"
15 #include "Thyra_DefaultSerialDenseLinearOpWithSolveFactory.hpp"
16 #include "Thyra_DetachedMultiVectorView.hpp"
17 #include "Thyra_DetachedVectorView.hpp"
18 #include "Thyra_MultiVectorStdOps.hpp"
19 #include "Thyra_VectorStdOps.hpp"
20 #include "Thyra_PreconditionerBase.hpp"
21 
22 // Epetra support
24 #include "Thyra_EpetraLinearOp.hpp"
25 #include "Thyra_get_Epetra_Operator.hpp"
26 #include "Epetra_Comm.h"
27 #include "Epetra_Map.h"
28 #include "Epetra_Vector.h"
29 #include "Epetra_Import.h"
30 #include "Epetra_CrsGraph.h"
31 #include "Epetra_CrsMatrix.h"
32 
33 namespace Tempus_Test {
34 
35 // Constructor
36 
37 template <class Scalar>
39  const int num_global_elements, const Scalar z_min,
40  const Scalar z_max, const Scalar a, const Scalar k)
41  : comm_(comm),
42  num_global_elements_(num_global_elements),
43  z_min_(z_min),
44  z_max_(z_max),
45  a_(a),
46  k_(k),
47  showGetInvalidArg_(false)
48 {
49  using Teuchos::RCP;
50  using Teuchos::rcp;
51  using ::Thyra::VectorBase;
52  typedef ::Thyra::ModelEvaluatorBase MEB;
54 
56 
57  const int num_nodes = num_global_elements_ + 1;
58 
59  // owned space
60  x_owned_map_ = rcp(new Epetra_Map(num_nodes, 0, *comm_));
62 
63  // ghosted space
64  if (comm_->NumProc() == 1) {
66  }
67  else {
68  int OverlapNumMyElements;
69  int OverlapMinMyGID;
70  OverlapNumMyElements = x_owned_map_->NumMyElements() + 2;
71  if ((comm_->MyPID() == 0) || (comm_->MyPID() == (comm_->NumProc() - 1)))
72  OverlapNumMyElements--;
73 
74  if (comm_->MyPID() == 0)
75  OverlapMinMyGID = x_owned_map_->MinMyGID();
76  else
77  OverlapMinMyGID = x_owned_map_->MinMyGID() - 1;
78 
79  int* OverlapMyGlobalElements = new int[OverlapNumMyElements];
80 
81  for (int i = 0; i < OverlapNumMyElements; i++)
82  OverlapMyGlobalElements[i] = OverlapMinMyGID + i;
83 
84  x_ghosted_map_ = Teuchos::rcp(new Epetra_Map(
85  -1, OverlapNumMyElements, OverlapMyGlobalElements, 0, *comm_));
86 
87  delete[] OverlapMyGlobalElements;
88  }
89 
90  importer_ = Teuchos::rcp(new Epetra_Import(*x_ghosted_map_, *x_owned_map_));
91 
92  // residual space
95 
96  x0_ = ::Thyra::createMember(x_space_);
97  V_S(x0_.ptr(), ST::zero());
98 
99  // set_p(Teuchos::tuple<Scalar>(p0, p1)());
100  // set_x0(Teuchos::tuple<Scalar>(x0, x1)());
101 
102  // Initialize the graph for W CrsMatrix object
103  W_graph_ = createGraph();
104 
105  // Create the nodal coordinates
106  {
107  node_coordinates_ = Teuchos::rcp(new Epetra_Vector(*x_owned_map_));
108  Scalar length = z_max_ - z_min_;
109  Scalar dx = length / ((double)num_global_elements_ - 1);
110  for (int i = 0; i < x_owned_map_->NumMyElements(); i++) {
111  (*node_coordinates_)[i] =
112  z_min_ + dx * ((double)x_owned_map_->MinMyGID() + i);
113  }
114  }
115 
116  MEB::InArgsSetup<Scalar> inArgs;
117  inArgs.setModelEvalDescription(this->description());
118  inArgs.setSupports(MEB::IN_ARG_t);
119  inArgs.setSupports(MEB::IN_ARG_x);
120  inArgs.setSupports(MEB::IN_ARG_beta);
121  inArgs.setSupports(MEB::IN_ARG_x_dot);
122  inArgs.setSupports(MEB::IN_ARG_alpha);
123  prototypeInArgs_ = inArgs;
124 
125  MEB::OutArgsSetup<Scalar> outArgs;
126  outArgs.setModelEvalDescription(this->description());
127  outArgs.setSupports(MEB::OUT_ARG_f);
128  outArgs.setSupports(MEB::OUT_ARG_W);
129  outArgs.setSupports(MEB::OUT_ARG_W_op);
130  outArgs.setSupports(MEB::OUT_ARG_W_prec);
131  // outArgs.set_W_properties(DerivativeProperties(
132  // DERIV_LINEARITY_NONCONST
133  // ,DERIV_RANK_FULL
134  // ,true // supportsAdjoint
135  // ));
136  prototypeOutArgs_ = outArgs;
137 
138  // Setup nominal values
139  nominalValues_ = inArgs;
140  nominalValues_.set_x(x0_);
141  auto x_dot_init = Thyra::createMember(this->get_x_space());
142  Thyra::put_scalar(Scalar(0.0), x_dot_init.ptr());
143  nominalValues_.set_x_dot(x_dot_init);
144 }
145 
146 // Initializers/Accessors
147 
148 template <class Scalar>
150 {
152 
153  // Create the shell for the
154  W_graph = Teuchos::rcp(new Epetra_CrsGraph(Copy, *x_owned_map_, 5));
155 
156  // Declare required variables
157  int row;
158  int column;
159  int OverlapNumMyElements = x_ghosted_map_->NumMyElements();
160 
161  // Loop Over # of Finite Elements on Processor
162  for (int ne = 0; ne < OverlapNumMyElements - 1; ne++) {
163  // Loop over Nodes in Element
164  for (int i = 0; i < 2; i++) {
165  row = x_ghosted_map_->GID(ne + i);
166 
167  // Loop over Trial Functions
168  for (int j = 0; j < 2; j++) {
169  // If this row is owned by current processor, add the index
170  if (x_owned_map_->MyGID(row)) {
171  column = x_ghosted_map_->GID(ne + j);
172  W_graph->InsertGlobalIndices(row, 1, &column);
173  }
174  }
175  }
176  }
177  W_graph->FillComplete();
178  return W_graph;
179 }
180 
181 template <class Scalar>
183 {
184 #ifdef TEUCHOS_DEBUG
185  TEUCHOS_ASSERT_EQUALITY(x_space_->dim(), x0_in.size());
186 #endif
188  x0.sv().values()().assign(x0_in);
189 }
190 
191 template <class Scalar>
192 void CDR_Model<Scalar>::setShowGetInvalidArgs(bool showGetInvalidArg)
193 {
194  showGetInvalidArg_ = showGetInvalidArg;
195 }
196 
197 template <class Scalar>
199  const Teuchos::RCP<const ::Thyra::LinearOpWithSolveFactoryBase<Scalar> >&
200  W_factory)
201 {
202  W_factory_ = W_factory;
203 }
204 
205 // Public functions overridden from ModelEvaluator
206 
207 template <class Scalar>
210 {
211  return x_space_;
212 }
213 
214 template <class Scalar>
217 {
218  return f_space_;
219 }
220 
221 template <class Scalar>
223  const
224 {
225  return nominalValues_;
226 }
227 
228 template <class Scalar>
231 {
233  this->get_W_factory();
234 
235  TEUCHOS_TEST_FOR_EXCEPTION(is_null(W_factory), std::runtime_error,
236  "W_factory in CDR_Model has a null W_factory!");
237 
238  Teuchos::RCP<Thyra::LinearOpBase<double> > matrix = this->create_W_op();
239 
241  Thyra::linearOpWithSolve<double>(*W_factory, matrix);
242 
243  return W;
244 }
245 
246 template <class Scalar>
248  const
249 {
251  Teuchos::rcp(new Epetra_CrsMatrix(::Copy, *W_graph_));
252 
253  return Thyra::nonconstEpetraLinearOp(W_epetra);
254 }
255 
256 template <class Scalar>
259 {
261  Teuchos::rcp(new Epetra_CrsMatrix(::Copy, *W_graph_));
262 
264  Thyra::nonconstEpetraLinearOp(W_epetra);
265 
268 
269  prec->initializeRight(W_op);
270 
271  return prec;
272 }
273 
274 template <class Scalar>
277 {
278  return W_factory_;
279 }
280 
281 template <class Scalar>
283  const
284 {
285  return prototypeInArgs_;
286 }
287 
288 // Private functions overridden from ModelEvaluatorDefaultBase
289 
290 template <class Scalar>
293 {
294  return prototypeOutArgs_;
295 }
296 
297 template <class Scalar>
300  const Thyra::ModelEvaluatorBase::OutArgs<Scalar>& outArgs) const
301 {
302  using Teuchos::RCP;
303  using Teuchos::rcp;
304  using Teuchos::rcp_dynamic_cast;
305 
306  TEUCHOS_ASSERT(nonnull(inArgs.get_x()));
307  TEUCHOS_ASSERT(nonnull(inArgs.get_x_dot()));
308 
309  // const Thyra::ConstDetachedVectorView<Scalar> x(inArgs.get_x());
310 
311  const RCP<Thyra::VectorBase<Scalar> > f_out = outArgs.get_f();
312  const RCP<Thyra::LinearOpBase<Scalar> > W_out = outArgs.get_W_op();
313  const RCP<Thyra::PreconditionerBase<Scalar> > W_prec_out =
314  outArgs.get_W_prec();
315 
316  if (nonnull(f_out) || nonnull(W_out) || nonnull(W_prec_out)) {
317  // ****************
318  // Get the underlying epetra objects
319  // ****************
320 
322  if (nonnull(f_out)) {
323  f = Thyra::get_Epetra_Vector(*f_owned_map_, outArgs.get_f());
324  }
325 
327  if (nonnull(W_out)) {
329  J = rcp_dynamic_cast<Epetra_CrsMatrix>(W_epetra);
331  }
332 
333  RCP<Epetra_CrsMatrix> M_inv;
334  if (nonnull(W_prec_out)) {
335  RCP<Epetra_Operator> M_epetra =
336  Thyra::get_Epetra_Operator(*(W_prec_out->getNonconstRightPrecOp()));
337  M_inv = rcp_dynamic_cast<Epetra_CrsMatrix>(M_epetra);
338  TEUCHOS_ASSERT(nonnull(M_inv));
339  J_diagonal_ = Teuchos::rcp(new Epetra_Vector(*x_owned_map_));
340  }
341 
342  // ****************
343  // Create ghosted objects
344  // ****************
345 
346  // Set the boundary condition directly. Works for both x and xDot solves.
347  if (comm_->MyPID() == 0) {
349  Teuchos::rcp_const_cast<Thyra::VectorBase<Scalar> >(inArgs.get_x());
350  (*Thyra::get_Epetra_Vector(*x_owned_map_, x))[0] = 1.0;
351  }
352 
353  if (is_null(u_ptr))
354  u_ptr = Teuchos::rcp(new Epetra_Vector(*x_ghosted_map_));
355 
356  u_ptr->Import(*(Thyra::get_Epetra_Vector(*x_owned_map_, inArgs.get_x())),
357  *importer_, Insert);
358 
359  if (is_null(u_dot_ptr))
360  u_dot_ptr = Teuchos::rcp(new Epetra_Vector(*x_ghosted_map_));
361 
362  u_dot_ptr->Import(
363  *(Thyra::get_Epetra_Vector(*x_owned_map_, inArgs.get_x_dot())),
364  *importer_, Insert);
365 
366  if (is_null(x_ptr)) {
367  x_ptr = Teuchos::rcp(new Epetra_Vector(*x_ghosted_map_));
368  x_ptr->Import(*node_coordinates_, *importer_, Insert);
369  }
370 
371  Epetra_Vector& u = *u_ptr;
372  Epetra_Vector& u_dot = *u_dot_ptr;
373  Epetra_Vector& x = *x_ptr;
374 
375  int ierr = 0;
376  int OverlapNumMyElements = x_ghosted_map_->NumMyElements();
377 
378  double xx[2];
379  double uu[2];
380  double uu_dot[2];
381  Basis basis;
382  const auto alpha = inArgs.get_alpha();
383  const auto beta = inArgs.get_beta();
384 
385  // Zero out the objects that will be filled
386  if (nonnull(f)) f->PutScalar(0.0);
387  if (nonnull(J)) J->PutScalar(0.0);
388  if (nonnull(M_inv)) M_inv->PutScalar(0.0);
389 
390  // Loop Over # of Finite Elements on Processor
391  for (int ne = 0; ne < OverlapNumMyElements - 1; ne++) {
392  // Loop Over Gauss Points
393  for (int gp = 0; gp < 2; gp++) {
394  // Get the solution and coordinates at the nodes
395  xx[0] = x[ne];
396  xx[1] = x[ne + 1];
397  uu[0] = u[ne];
398  uu[1] = u[ne + 1];
399  uu_dot[0] = u_dot[ne];
400  uu_dot[1] = u_dot[ne + 1];
401  // Calculate the basis function at the gauss point
402  basis.computeBasis(gp, xx, uu, uu_dot);
403 
404  // Loop over Nodes in Element
405  for (int i = 0; i < 2; i++) {
406  int row = x_ghosted_map_->GID(ne + i);
407  // printf("Proc=%d GlobalRow=%d LocalRow=%d Owned=%d\n",
408  // MyPID, row, ne+i,x_owned_map_.MyGID(row));
409  if (x_owned_map_->MyGID(row)) {
410  if (nonnull(f)) {
411  (*f)[x_owned_map_->LID(x_ghosted_map_->GID(ne + i))] +=
412  +basis.wt * basis.dz *
413  (basis.uu_dot * basis.phi[i] // transient
414  + (a_ / basis.dz * basis.duu * basis.phi[i] // convection
415  + 1.0 / (basis.dz * basis.dz)) *
416  basis.duu * basis.dphide[i] // diffusion
417  + k_ * basis.uu * basis.uu * basis.phi[i]); // source
418  }
419  }
420  // Loop over Trial Functions
421  if (nonnull(J)) {
422  for (int j = 0; j < 2; j++) {
423  if (x_owned_map_->MyGID(row)) {
424  int column = x_ghosted_map_->GID(ne + j);
425  double jac = basis.wt * basis.dz *
426  (alpha * basis.phi[i] * basis.phi[j] // transient
427  + beta * (+a_ / basis.dz * basis.dphide[j] *
428  basis.phi[i] // convection
429  + (1.0 / (basis.dz * basis.dz)) *
430  basis.dphide[j] *
431  basis.dphide[i] // diffusion
432  + 2.0 * k_ * basis.uu * basis.phi[j] *
433  basis.phi[i] // source
434  ));
435  ierr = J->SumIntoGlobalValues(row, 1, &jac, &column);
436  }
437  }
438  }
439  if (nonnull(M_inv)) {
440  for (int j = 0; j < 2; j++) {
441  if (x_owned_map_->MyGID(row)) {
442  int column = x_ghosted_map_->GID(ne + j);
443  // The prec will be the diagonal of J. No need to assemble the
444  // other entries
445  if (row == column) {
446  double jac =
447  basis.wt * basis.dz *
448  (alpha * basis.phi[i] * basis.phi[j] // transient
449  + beta * (+a_ / basis.dz * basis.dphide[j] *
450  basis.phi[i] // convection
451  + (1.0 / (basis.dz * basis.dz)) *
452  basis.dphide[j] *
453  basis.dphide[i] // diffusion
454  + 2.0 * k_ * basis.uu * basis.phi[j] *
455  basis.phi[i] // source
456  ));
457  ierr = M_inv->SumIntoGlobalValues(row, 1, &jac, &column);
458  }
459  }
460  }
461  }
462  }
463  }
464  }
465 
466  // Insert Boundary Conditions and modify Jacobian and function (F)
467  // U(0)=1
468  if (comm_->MyPID() == 0) {
469  if (nonnull(f))
470  (*f)[0] = 0.0; // Setting BC above and zero residual here works for x
471  // and xDot solves.
472  //(*f)[0]= u[0] - 1.0; // BC equation works for x solves.
473  if (nonnull(J)) {
474  int column = 0;
475  double jac = 1.0;
476  ierr = J->ReplaceGlobalValues(0, 1, &jac, &column);
477  column = 1;
478  jac = 0.0;
479  ierr = J->ReplaceGlobalValues(0, 1, &jac, &column);
480  }
481  if (nonnull(M_inv)) {
482  int column = 0;
483  double jac = 1.0;
484  ierr = M_inv->ReplaceGlobalValues(0, 1, &jac, &column);
485  column = 1;
486  jac = 0.0;
487  ierr = M_inv->ReplaceGlobalValues(0, 1, &jac, &column);
488  }
489  }
490 
491  if (nonnull(J)) J->FillComplete();
492 
493  if (nonnull(M_inv)) {
494  // invert the Jacobian diagonal for the preconditioner
495  M_inv->ExtractDiagonalCopy(*J_diagonal_);
496 
497  for (int i = 0; i < J_diagonal_->MyLength(); ++i)
498  (*J_diagonal_)[i] = 1.0 / ((*J_diagonal_)[i]);
499 
500  M_inv->PutScalar(0.0);
501  M_inv->ReplaceDiagonalValues(*J_diagonal_);
502  M_inv->FillComplete();
503  }
504 
505  TEUCHOS_ASSERT(ierr > -1);
506  }
507 }
508 
509 //====================================================================
510 // Basis vector
511 
512 // Constructor
514  : uu(0.0),
515  zz(0.0),
516  duu(0.0),
517  eta(0.0),
518  wt(0.0),
519  dz(0.0),
520  uu_dot(0.0),
521  duu_dot(0.0)
522 {
523  phi = new double[2];
524  dphide = new double[2];
525 }
526 
527 // Destructor
529 {
530  delete[] phi;
531  delete[] dphide;
532 }
533 
534 // Calculates a linear 1D basis
535 void Basis::computeBasis(int gp, double* z, double* u, double* u_dot)
536 {
537  int N = 2;
538  if (gp == 0) {
539  eta = -1.0 / sqrt(3.0);
540  wt = 1.0;
541  }
542  if (gp == 1) {
543  eta = 1.0 / sqrt(3.0);
544  wt = 1.0;
545  }
546 
547  // Calculate basis function and derivatives at nodel pts
548  phi[0] = (1.0 - eta) / 2.0;
549  phi[1] = (1.0 + eta) / 2.0;
550  dphide[0] = -0.5;
551  dphide[1] = 0.5;
552 
553  // Caculate basis function and derivative at GP.
554  dz = 0.5 * (z[1] - z[0]);
555  zz = 0.0;
556  uu = 0.0;
557  duu = 0.0;
558  uu_dot = 0.0;
559  duu_dot = 0.0;
560  for (int i = 0; i < N; i++) {
561  zz += z[i] * phi[i];
562  uu += u[i] * phi[i];
563  duu += u[i] * dphide[i];
564  if (u_dot) {
565  uu_dot += u_dot[i] * phi[i];
566  duu_dot += u_dot[i] * dphide[i];
567  }
568  }
569 
570  return;
571 }
572 
573 } // namespace Tempus_Test
574 
575 #endif
::Thyra::ModelEvaluatorBase::InArgs< Scalar > createInArgs() const
void evalModelImpl(const ::Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const ::Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
Teuchos::RCP< const Epetra_Map > x_ghosted_map_
bool is_null(const boost::shared_ptr< T > &p)
Teuchos::RCP< const ::Thyra::LinearOpWithSolveFactoryBase< Scalar > > get_W_factory() const
CDR_Model(const Teuchos::RCP< const Epetra_Comm > &comm, const int num_global_elements, const Scalar z_min, const Scalar z_max, const Scalar a, const Scalar k)
::Thyra::ModelEvaluatorBase::OutArgs< Scalar > prototypeOutArgs_
virtual Teuchos::RCP< Epetra_CrsGraph > createGraph()
RCP< const VectorBase< Scalar > > get_x_dot() const
RCP< const VectorSpaceBase< double > > create_VectorSpace(const RCP< const Epetra_Map > &epetra_map)
Teuchos::RCP< Epetra_CrsGraph > W_graph_
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< Epetra_Vector > node_coordinates_
Teuchos::RCP< const Epetra_Map > f_owned_map_
size_type size() const
Evaluation< VectorBase< Scalar > > get_f() const
Teuchos::RCP< ::Thyra::VectorBase< Scalar > > x0_
Teuchos::RCP< ::Thyra::PreconditionerBase< Scalar > > create_W_prec() const
void setShowGetInvalidArgs(bool showGetInvalidArg)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< const Epetra_Map > x_owned_map_
Ptr< T > ptr() const
::Thyra::ModelEvaluatorBase::InArgs< Scalar > nominalValues_
void set_W_factory(const Teuchos::RCP< const ::Thyra::LinearOpWithSolveFactoryBase< Scalar > > &W_factory)
void set_x0(const Teuchos::ArrayView< const Scalar > &x0)
RCP< Epetra_Vector > get_Epetra_Vector(const Epetra_Map &map, const RCP< VectorBase< double > > &v)
virtual std::string description() const
void computeBasis(int gp, double *z, double *u, double *u_dot=nullptr)
Teuchos::RCP< const ::Thyra::VectorSpaceBase< Scalar > > x_space_
Teuchos::RCP< ::Thyra::LinearOpBase< Scalar > > create_W_op() const
Teuchos::RCP< Thyra::LinearOpWithSolveBase< double > > create_W() const
::Thyra::ModelEvaluatorBase::OutArgs< Scalar > createOutArgsImpl() const
bool nonnull(const boost::shared_ptr< T > &p)
Teuchos::RCP< const ::Thyra::VectorSpaceBase< Scalar > > get_f_space() const
Teuchos::RCP< Epetra_Operator > get_Epetra_Operator(LinearOpBase< Scalar > &)
::Thyra::ModelEvaluatorBase::InArgs< Scalar > prototypeInArgs_
RCP< PreconditionerBase< Scalar > > get_W_prec() const
Teuchos::RCP< const ::Thyra::VectorSpaceBase< Scalar > > get_x_space() const
#define TEUCHOS_ASSERT(assertion_test)
#define TEUCHOS_ASSERT_EQUALITY(val1, val2)
const Teuchos::RCP< const Epetra_Comm > comm_
Teuchos::RCP< const ::Thyra::VectorSpaceBase< Scalar > > f_space_
RCP< LinearOpBase< Scalar > > get_W_op() const
RCP< const VectorBase< Scalar > > get_x() const
::Thyra::ModelEvaluatorBase::InArgs< Scalar > getNominalValues() const
Teuchos::RCP< const Epetra_Import > importer_