10 #ifndef NOX_THYRA_MODEL_EVALUATOR_1DFEM_DEF_HPP
11 #define NOX_THYRA_MODEL_EVALUATOR_1DFEM_DEF_HPP
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"
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"
33 namespace Tempus_Test {
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)
42 num_global_elements_(num_global_elements),
47 showGetInvalidArg_(false)
51 using ::Thyra::VectorBase;
52 typedef ::Thyra::ModelEvaluatorBase MEB;
64 if (
comm_->NumProc() == 1) {
68 int OverlapNumMyElements;
70 OverlapNumMyElements =
x_owned_map_->NumMyElements() + 2;
71 if ((
comm_->MyPID() == 0) || (
comm_->MyPID() == (
comm_->NumProc() - 1)))
72 OverlapNumMyElements--;
74 if (
comm_->MyPID() == 0)
79 int* OverlapMyGlobalElements =
new int[OverlapNumMyElements];
81 for (
int i = 0; i < OverlapNumMyElements; i++)
82 OverlapMyGlobalElements[i] = OverlapMinMyGID + i;
85 -1, OverlapNumMyElements, OverlapMyGlobalElements, 0, *
comm_));
87 delete[] OverlapMyGlobalElements;
97 V_S(
x0_.
ptr(), ST::zero());
110 for (
int i = 0; i <
x_owned_map_->NumMyElements(); i++) {
111 (*node_coordinates_)[i] =
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);
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);
141 auto x_dot_init = Thyra::createMember(this->
get_x_space());
142 Thyra::put_scalar(Scalar(0.0), x_dot_init.ptr());
148 template <
class Scalar>
159 int OverlapNumMyElements = x_ghosted_map_->NumMyElements();
162 for (
int ne = 0; ne < OverlapNumMyElements - 1; ne++) {
164 for (
int i = 0; i < 2; i++) {
165 row = x_ghosted_map_->GID(ne + i);
168 for (
int j = 0; j < 2; j++) {
170 if (x_owned_map_->MyGID(row)) {
171 column = x_ghosted_map_->GID(ne + j);
172 W_graph->InsertGlobalIndices(row, 1, &column);
177 W_graph->FillComplete();
181 template <
class Scalar>
188 x0.sv().values()().assign(x0_in);
191 template <
class Scalar>
194 showGetInvalidArg_ = showGetInvalidArg;
197 template <
class Scalar>
199 const Teuchos::RCP<const ::Thyra::LinearOpWithSolveFactoryBase<Scalar> >&
202 W_factory_ = W_factory;
207 template <
class Scalar>
214 template <
class Scalar>
221 template <
class Scalar>
225 return nominalValues_;
228 template <
class Scalar>
233 this->get_W_factory();
236 "W_factory in CDR_Model has a null W_factory!");
241 Thyra::linearOpWithSolve<double>(*W_factory, matrix);
246 template <
class Scalar>
253 return Thyra::nonconstEpetraLinearOp(W_epetra);
256 template <
class Scalar>
264 Thyra::nonconstEpetraLinearOp(W_epetra);
269 prec->initializeRight(W_op);
274 template <
class Scalar>
281 template <
class Scalar>
285 return prototypeInArgs_;
290 template <
class Scalar>
294 return prototypeOutArgs_;
297 template <
class Scalar>
304 using Teuchos::rcp_dynamic_cast;
329 J = rcp_dynamic_cast<Epetra_CrsMatrix>(W_epetra);
337 M_inv = rcp_dynamic_cast<Epetra_CrsMatrix>(M_epetra);
339 J_diagonal_ =
Teuchos::rcp(
new Epetra_Vector(*x_owned_map_));
347 if (comm_->MyPID() == 0) {
354 u_ptr =
Teuchos::rcp(
new Epetra_Vector(*x_ghosted_map_));
360 u_dot_ptr =
Teuchos::rcp(
new Epetra_Vector(*x_ghosted_map_));
367 x_ptr =
Teuchos::rcp(
new Epetra_Vector(*x_ghosted_map_));
368 x_ptr->Import(*node_coordinates_, *importer_, Insert);
371 Epetra_Vector& u = *u_ptr;
372 Epetra_Vector& u_dot = *u_dot_ptr;
373 Epetra_Vector& x = *x_ptr;
376 int OverlapNumMyElements = x_ghosted_map_->NumMyElements();
383 const auto beta = inArgs.
get_beta();
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);
391 for (
int ne = 0; ne < OverlapNumMyElements - 1; ne++) {
393 for (
int gp = 0; gp < 2; gp++) {
399 uu_dot[0] = u_dot[ne];
400 uu_dot[1] = u_dot[ne + 1];
405 for (
int i = 0; i < 2; i++) {
406 int row = x_ghosted_map_->GID(ne + i);
409 if (x_owned_map_->MyGID(row)) {
411 (*f)[x_owned_map_->LID(x_ghosted_map_->GID(ne + i))] +=
412 +basis.
wt * basis.
dz *
414 + (a_ / basis.
dz * basis.
duu * basis.
phi[i]
415 + 1.0 / (basis.
dz * basis.
dz)) *
417 + k_ * basis.
uu * basis.
uu * basis.
phi[i]);
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]
427 + beta * (+a_ / basis.
dz * basis.
dphide[j] *
429 + (1.0 / (basis.
dz * basis.
dz)) *
432 + 2.0 * k_ * basis.
uu * basis.
phi[j] *
435 ierr = J->SumIntoGlobalValues(row, 1, &jac, &column);
440 for (
int j = 0; j < 2; j++) {
441 if (x_owned_map_->MyGID(row)) {
442 int column = x_ghosted_map_->GID(ne + j);
447 basis.
wt * basis.
dz *
448 (alpha * basis.
phi[i] * basis.
phi[j]
449 + beta * (+a_ / basis.
dz * basis.
dphide[j] *
451 + (1.0 / (basis.
dz * basis.
dz)) *
454 + 2.0 * k_ * basis.
uu * basis.
phi[j] *
457 ierr = M_inv->SumIntoGlobalValues(row, 1, &jac, &column);
468 if (comm_->MyPID() == 0) {
476 ierr = J->ReplaceGlobalValues(0, 1, &jac, &column);
479 ierr = J->ReplaceGlobalValues(0, 1, &jac, &column);
484 ierr = M_inv->ReplaceGlobalValues(0, 1, &jac, &column);
487 ierr = M_inv->ReplaceGlobalValues(0, 1, &jac, &column);
491 if (
nonnull(J)) J->FillComplete();
495 M_inv->ExtractDiagonalCopy(*J_diagonal_);
497 for (
int i = 0; i < J_diagonal_->MyLength(); ++i)
498 (*J_diagonal_)[i] = 1.0 / ((*J_diagonal_)[i]);
500 M_inv->PutScalar(0.0);
501 M_inv->ReplaceDiagonalValues(*J_diagonal_);
502 M_inv->FillComplete();
539 eta = -1.0 / sqrt(3.0);
543 eta = 1.0 / sqrt(3.0);
548 phi[0] = (1.0 -
eta) / 2.0;
549 phi[1] = (1.0 +
eta) / 2.0;
554 dz = 0.5 * (z[1] - z[0]);
560 for (
int i = 0; i < N; i++) {
565 uu_dot += u_dot[i] * phi[i];
566 duu_dot += u_dot[i] * dphide[i];
::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_
Evaluation< VectorBase< Scalar > > get_f() const
const int num_global_elements_
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_
::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_