Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
CDR_Model_Functors.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 TEMPUS_CDR_MODEL_FUNCTORS_HPP
11 #define TEMPUS_CDR_MODEL_FUNCTORS_HPP
12 
13 #include "Kokkos_Core.hpp"
14 #include "Tpetra_Details_OrdinalTraits.hpp"
15 
16 namespace Tempus_Test {
17 
18 template <class TpetraVectorType>
19 struct CoordFiller {
20  using SC = typename TpetraVectorType::impl_scalar_type;
21  using LO = typename TpetraVectorType::local_ordinal_type;
22  using GO = typename TpetraVectorType::global_ordinal_type;
23  using Map = typename TpetraVectorType::map_type;
24  using LocalMap = typename Map::local_map_type;
25  using DV = typename TpetraVectorType::dual_view_type;
26  using View = typename DV::t_dev;
27 
29  const SC zMin_;
30  const SC dx_;
31  const SC minGI_;
32 
33  CoordFiller(TpetraVectorType &coords, const SC zMin, const SC dx,
34  const GO minGI)
35  : coordsView_(coords.getLocalViewDevice(Tpetra::Access::ReadWrite)),
36  zMin_(zMin),
37  dx_(dx),
38  minGI_(minGI)
39  {
40  }
41 
42  KOKKOS_INLINE_FUNCTION
43  void operator()(const LO i) const
44  {
45  coordsView_(i, 0) = zMin_ + dx_ * static_cast<SC>(minGI_ + i);
46  }
47 };
48 // Finite Element Basis Object
49 template <class Scalar, class LO>
50 class TBasis {
51  public:
52  // Calculates the values of z and u at the specified Gauss point
53  KOKKOS_INLINE_FUNCTION
54  void computeBasis(LO gp, Scalar *z, Scalar *u, Scalar *u_dot)
55  {
56  if (gp == 0) {
57  eta = -1.0 / sqrt(3.0);
58  wt = 1.0;
59  }
60  if (gp == 1) {
61  eta = 1.0 / sqrt(3.0);
62  wt = 1.0;
63  }
64 
65  // Calculate basis function and derivatives at Gauss point
66  phi[0] = (1.0 - eta) / 2.0;
67  phi[1] = (1.0 + eta) / 2.0;
68  dphide[0] = -0.5;
69  dphide[1] = 0.5;
70 
71  // Caculate function and derivative approximations at GP.
72  dz = 0.5 * (z[1] - z[0]);
73  zz = 0.0;
74  uu = 0.0;
75  duu = 0.0;
76  uu_dot = 0.0;
77  duu_dot = 0.0;
78 
79  for (LO i = 0; i < 2; i++) {
80  zz += z[i] * phi[i];
81  uu += u[i] * phi[i];
82  duu += u[i] * dphide[i];
83  uu_dot += u_dot[i] * phi[i];
84  duu_dot += u_dot[i] * dphide[i];
85  }
86  }
87 
88  public:
89  // Variables that are calculated at the Gauss point
90  Scalar phi[2];
91  Scalar dphide[2];
92  Scalar uu = 0.0;
93  Scalar zz = 0.0;
94  Scalar duu = 0.0;
95  Scalar eta = 0.0;
96  Scalar wt = 0.0;
97  Scalar dz = 0.0;
98 
99  // These are only needed for transient
100  Scalar uu_dot = 0.0;
101  Scalar duu_dot = 0.0;
102 };
103 
104 //==================================================================
105 
106 template <class TpetraVectorType, class TpetraMatrixType>
108  using SC = typename TpetraVectorType::impl_scalar_type;
109  using LO = typename TpetraVectorType::local_ordinal_type;
110  using Map = typename TpetraVectorType::map_type;
111  using LocalMap = typename Map::local_map_type;
112  using DV = typename TpetraVectorType::dual_view_type;
113  using ConstView = typename DV::t_dev::const_type;
114  using LocalMat = typename TpetraMatrixType::local_matrix_device_type;
115 
122  const int myRank_;
123  const SC a_;
124  const SC k_;
125  const SC alpha_;
126  const SC beta_;
127 
128  JacobianEvaluatorFunctor(const TpetraMatrixType &J, const TpetraVectorType &x,
129  const TpetraVectorType &u,
130  const TpetraVectorType &uDot, const int &myRank,
131  const SC a, const SC k, const SC alpha,
132  const SC beta)
133  : jLocal_(J.getLocalMatrixDevice()),
134  xView_(x.getLocalViewDevice(Tpetra::Access::ReadOnly)),
135  uView_(u.getLocalViewDevice(Tpetra::Access::ReadOnly)),
136  uDotView_(uDot.getLocalViewDevice(Tpetra::Access::ReadOnly)),
137  rowMap_(J.getRowMap()->getLocalMap()),
138  colMap_(J.getColMap()->getLocalMap()),
139  myRank_(myRank),
140  a_(a),
141  k_(k),
142  alpha_(alpha),
143  beta_(beta)
144  {
145  }
146 
147  // Adds the contribution from element ne to the Jacobian matrix
148  KOKKOS_INLINE_FUNCTION
149  void operator()(const LO ne) const
150  {
151  const auto invalid = Tpetra::Details::OrdinalTraits<LO>::invalid();
152  // Get the solution and coordinates at the nodes
153  SC xx[2];
154  xx[0] = xView_(ne, 0);
155  xx[1] = xView_(ne + 1, 0);
156 
157  SC uu[2];
158  uu[0] = uView_(ne, 0);
159  uu[1] = uView_(ne + 1, 0);
160 
161  SC uuDot[2];
162  uuDot[0] = uDotView_(ne, 0);
163  uuDot[1] = uDotView_(ne + 1, 0);
164 
165  TBasis<SC, LO> basis;
166  // Loop Over Gauss Points
167  for (LO gp = 0; gp < 2; ++gp) {
168  // Calculate the basis function at the gauss point
169  basis.computeBasis(gp, xx, uu, uuDot);
170 
171  // Loop over nodes in element
172  for (LO i = 0; i < 2; ++i) {
173  auto localRow =
174  rowMap_.getLocalElement(colMap_.getGlobalElement(ne + i));
175 
176  if (localRow != invalid) {
177  // Loop over trial functions
178  for (LO j = 0; j < 2; ++j) {
179  const auto localColumn = ne + j;
180  double jac =
181  basis.wt * basis.dz *
182  (alpha_ * basis.phi[i] * basis.phi[j] // transient
183  + beta_ * (+a_ / basis.dz * basis.dphide[j] *
184  basis.phi[i] // convection
185  + (1.0 / (basis.dz * basis.dz)) * basis.dphide[j] *
186  basis.dphide[i] // diffusion
187  + 2.0 * k_ * basis.uu * basis.phi[j] *
188  basis.phi[i] // source
189  ));
190 
191  jLocal_.sumIntoValues(localRow, &localColumn, 1, &jac, false, true);
192  }
193  }
194  }
195 
196  // Correct for Dirichlet BCs
197  if ((myRank_ == 0) && (ne == 0)) {
198  LO row = 0;
199  LO col = 0;
200  SC val = 1.0;
201  jLocal_.replaceValues(row, &col, 1, &val);
202 
203  col = 1;
204  val = 0.0;
205  jLocal_.replaceValues(row, &col, 1, &val);
206  }
207  }
208  }
209 };
210 
211 //==================================================================
212 
213 template <class TpetraVectorType, class TpetraMatrixType>
215  using SC = typename TpetraVectorType::impl_scalar_type;
216  using LO = typename TpetraVectorType::local_ordinal_type;
217  using Map = typename TpetraVectorType::map_type;
218  using LocalMap = typename Map::local_map_type;
219  using DV = typename TpetraVectorType::dual_view_type;
220  using ConstView = typename DV::t_dev::const_type;
221  using LocalMat = typename TpetraMatrixType::local_matrix_device_type;
222 
229  const int myRank_;
230  const SC a_;
231  const SC k_;
232  const SC alpha_;
233  const SC beta_;
234 
235  PreconditionerEvaluatorFunctor(const TpetraMatrixType &M,
236  const TpetraVectorType &x,
237  const TpetraVectorType &u,
238  const TpetraVectorType &uDot,
239  const int &myRank, const SC a, const SC k,
240  const SC alpha, const SC beta)
241  : mLocal_(M.getLocalMatrixDevice()),
242  xView_(x.getLocalViewDevice(Tpetra::Access::ReadOnly)),
243  uView_(u.getLocalViewDevice(Tpetra::Access::ReadOnly)),
244  uDotView_(uDot.getLocalViewDevice(Tpetra::Access::ReadOnly)),
245  rowMap_(M.getRowMap()->getLocalMap()),
246  colMap_(M.getColMap()->getLocalMap()),
247  myRank_(myRank),
248  a_(a),
249  k_(k),
250  alpha_(alpha),
251  beta_(beta)
252  {
253  }
254 
255  // Adds the contribution from element ne to the preconditioner matrix
256  KOKKOS_INLINE_FUNCTION
257  void operator()(const LO ne) const
258  {
259  const auto invalid = Tpetra::Details::OrdinalTraits<LO>::invalid();
260  // Get the solution and coordinates at the nodes
261  SC xx[2];
262  xx[0] = xView_(ne, 0);
263  xx[1] = xView_(ne + 1, 0);
264 
265  SC uu[2];
266  uu[0] = uView_(ne, 0);
267  uu[1] = uView_(ne + 1, 0);
268 
269  SC uuDot[2];
270  uuDot[0] = uDotView_(ne, 0);
271  uuDot[1] = uDotView_(ne + 1, 0);
272 
273  TBasis<SC, LO> basis;
274 
275  // Loop Over Gauss Points
276  for (LO gp = 0; gp < 2; ++gp) {
277  // Calculate the basis function at the gauss point
278  basis.computeBasis(gp, xx, uu, uuDot);
279 
280  // Loop over nodes in element
281  for (LO i = 0; i < 2; ++i) {
282  auto localRow =
283  rowMap_.getLocalElement(colMap_.getGlobalElement(ne + i));
284 
285  // Loop over trial functions
286  if (localRow != invalid) {
287  for (LO j = 0; j < 2; ++j) {
288  const auto localColumn = ne + j;
289  if (rowMap_.getGlobalElement(localRow) ==
290  colMap_.getGlobalElement(localColumn)) {
291  auto value = basis.wt * basis.dz *
292  (alpha_ * basis.phi[i] * basis.phi[j] // transient
293  + beta_ * (+a_ / basis.dz * basis.dphide[j] *
294  basis.phi[i] // convection
295  + (1.0 / (basis.dz * basis.dz)) *
296  basis.dphide[j] *
297  basis.dphide[i] // diffusion
298  + 2.0 * k_ * basis.uu * basis.phi[j] *
299  basis.phi[i] // source
300  ));
301  mLocal_.sumIntoValues(localRow, &localColumn, 1, &value);
302  }
303  }
304  }
305  }
306  }
307 
308  // Correct for Dirichlet BCs
309  if ((myRank_ == 0) && (ne == 0)) {
310  LO row = 0;
311  LO column = 0;
312  SC value = 1.0;
313  mLocal_.replaceValues(row, &column, 1, &value);
314  }
315  }
316 };
317 
318 //==================================================================
319 
320 template <class TpetraVectorType>
322  using SC = typename TpetraVectorType::impl_scalar_type;
323  using LO = typename TpetraVectorType::local_ordinal_type;
324  using Map = typename TpetraVectorType::map_type;
325  using LocalMap = typename Map::local_map_type;
326  using DV = typename TpetraVectorType::dual_view_type;
327  using View = typename DV::t_dev;
328  using ConstView = typename DV::t_dev::const_type;
329 
330  const View fView_;
336  const int myRank_;
337  const SC a_;
338  const SC k_;
339 
340  DfDp2EvaluatorFunctor(TpetraVectorType &f, const TpetraVectorType &x,
341  const TpetraVectorType &u, const TpetraVectorType &uDot,
342  const int myRank, SC a, SC k)
343  : fView_(f.getLocalViewDevice(Tpetra::Access::ReadWrite)),
344  xView_(x.getLocalViewDevice(Tpetra::Access::ReadOnly)),
345  uView_(u.getLocalViewDevice(Tpetra::Access::ReadOnly)),
346  uDotView_(uDot.getLocalViewDevice(Tpetra::Access::ReadOnly)),
347  rowMap_(f.getMap()->getLocalMap()),
348  colMap_(u.getMap()->getLocalMap()),
349  myRank_(myRank),
350  a_(a),
351  k_(k)
352  {
353  }
354 
355  // Adds the contribution from element ne to the residual vector
356  KOKKOS_INLINE_FUNCTION
357  void operator()(const LO ne) const
358  {
359  const auto invalid = Tpetra::Details::OrdinalTraits<LO>::invalid();
360  // Get the solution and coordinates at the nodes
361  SC xx[2];
362  xx[0] = xView_(ne, 0);
363  xx[1] = xView_(ne + 1, 0);
364 
365  SC uu[2];
366  uu[0] = uView_(ne, 0);
367  uu[1] = uView_(ne + 1, 0);
368 
369  SC uuDot[2];
370  uuDot[0] = uDotView_(ne, 0);
371  uuDot[1] = uDotView_(ne + 1, 0);
372 
373  TBasis<SC, LO> basis;
374  // Loop Over Gauss Points
375  for (LO gp = 0; gp < 2; ++gp) {
376  // Calculate the basis function at the gauss point
377  basis.computeBasis(gp, xx, uu, uuDot);
378 
379  // Loop over nodes in element
380  for (LO i = 0; i < 2; ++i) {
381  auto localRow =
382  rowMap_.getLocalElement(colMap_.getGlobalElement(ne + i));
383  if (localRow != invalid) {
384  auto value =
385  basis.wt * basis.dz *
386  (basis.uu_dot * basis.phi[i] // transient
387  + (a_ / basis.dz * basis.duu * basis.phi[i] // convection
388  + 1.0 / (basis.dz * basis.dz)) *
389  basis.duu * basis.dphide[i] // diffusion
390  + k_ * basis.uu * basis.uu * basis.phi[i]); // source
391  Kokkos::atomic_add(&fView_(localRow, 0), value);
392  }
393  }
394  }
395 
396  // Correct for Dirichlet BCs
397  if ((myRank_ == 0) && (ne == 0)) {
398  SC value = 0.0;
399  Kokkos::atomic_exchange(&fView_(0, 0), value);
400  }
401  }
402 };
403 
404 } // namespace Tempus_Test
405 
406 #endif // TEMPUS_CDR_MODEL_FUNCTORS_HPP
typename TpetraVectorType::map_type Map
KOKKOS_INLINE_FUNCTION void operator()(const LO ne) const
typename Map::local_map_type LocalMap
typename Map::local_map_type LocalMap
typename TpetraVectorType::global_ordinal_type GO
typename TpetraVectorType::dual_view_type DV
typename DV::t_dev::const_type ConstView
typename DV::t_dev::const_type ConstView
KOKKOS_INLINE_FUNCTION void operator()(const LO i) const
typename TpetraVectorType::dual_view_type DV
typename TpetraVectorType::local_ordinal_type LO
DfDp2EvaluatorFunctor(TpetraVectorType &f, const TpetraVectorType &x, const TpetraVectorType &u, const TpetraVectorType &uDot, const int myRank, SC a, SC k)
JacobianEvaluatorFunctor(const TpetraMatrixType &J, const TpetraVectorType &x, const TpetraVectorType &u, const TpetraVectorType &uDot, const int &myRank, const SC a, const SC k, const SC alpha, const SC beta)
typename TpetraVectorType::dual_view_type DV
typename TpetraMatrixType::local_matrix_device_type LocalMat
typename TpetraVectorType::impl_scalar_type SC
typename TpetraVectorType::impl_scalar_type SC
KOKKOS_INLINE_FUNCTION void operator()(const LO ne) const
typename TpetraVectorType::local_ordinal_type LO
typename TpetraVectorType::impl_scalar_type SC
typename Map::local_map_type LocalMap
typename TpetraVectorType::map_type Map
PreconditionerEvaluatorFunctor(const TpetraMatrixType &M, const TpetraVectorType &x, const TpetraVectorType &u, const TpetraVectorType &uDot, const int &myRank, const SC a, const SC k, const SC alpha, const SC beta)
typename TpetraVectorType::map_type Map
CoordFiller(TpetraVectorType &coords, const SC zMin, const SC dx, const GO minGI)
typename TpetraVectorType::local_ordinal_type LO
typename TpetraVectorType::local_ordinal_type LO
typename TpetraVectorType::dual_view_type DV
KOKKOS_INLINE_FUNCTION void operator()(const LO ne) const
typename TpetraMatrixType::local_matrix_device_type LocalMat
typename TpetraVectorType::map_type Map
typename TpetraVectorType::impl_scalar_type SC
KOKKOS_INLINE_FUNCTION void computeBasis(LO gp, Scalar *z, Scalar *u, Scalar *u_dot)