Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_KLMatrixFreeOperator.cpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
11 #include "EpetraExt_BlockMultiVector.h"
12 #include "EpetraExt_BlockUtility.h"
13 
17  const Teuchos::RCP<const Stokhos::OrthogPolyBasis<int,double> >& sg_basis_,
19  const Teuchos::RCP<const Epetra_Map>& domain_base_map_,
20  const Teuchos::RCP<const Epetra_Map>& range_base_map_,
21  const Teuchos::RCP<const Epetra_Map>& domain_sg_map_,
22  const Teuchos::RCP<const Epetra_Map>& range_sg_map_,
24  label("Stokhos KL MatrixFree Operator"),
25  sg_comm(sg_comm_),
26  sg_basis(sg_basis_),
27  epetraCijk(epetraCijk_),
28  domain_base_map(domain_base_map_),
29  range_base_map(range_base_map_),
30  domain_sg_map(domain_sg_map_),
31  range_sg_map(range_sg_map_),
32  is_stoch_parallel(epetraCijk->isStochasticParallel()),
33  global_col_map(),
34  global_col_map_trans(),
35  stoch_col_map(epetraCijk->getStochasticColMap()),
36  col_importer(),
37  col_importer_trans(),
38  Cijk(epetraCijk->getParallelCijk()),
39  block_ops(),
40  scale_op(true),
41  include_mean(true),
42  useTranspose(false),
43  expansion_size(sg_basis->size()),
44  num_blocks(0),
45  input_col(),
46  input_col_trans(),
47  input_block(),
48  result_block(),
49  tmp(),
50  tmp_trans(),
51  k_begin(Cijk->k_begin()),
52  k_end(Cijk->k_end())
53 {
54  scale_op = params->get("Scale Operator by Inverse Basis Norms", true);
55  include_mean = params->get("Include Mean", true);
56 
57  // Compute maximum number of mat-vec's needed
58  if (!include_mean && index(k_begin) == 0)
59  ++k_begin;
60  int dim = sg_basis->dimension();
61  k_end = Cijk->find_k(dim+1);
62  max_num_mat_vec = 0;
63  for (Cijk_type::k_iterator k=k_begin; k!=k_end; ++k) {
64  int nj = Cijk->num_j(k);
65  if (max_num_mat_vec < nj)
66  max_num_mat_vec = nj;
67  }
68 
69  // Build up column map of SG operator
70  int num_col_blocks = expansion_size;
71  int num_row_blocks = expansion_size;
72  if (is_stoch_parallel) {
73 
74  // Build column map from base domain map. This will communicate
75  // stochastic components to column map, but not deterministic. It would
76  // be more efficient to do both, but the Epetra_Operator interface
77  // doesn't have the concept of a column map.
79  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(*domain_base_map,
81  *sg_comm));
83  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(*range_base_map,
85  *sg_comm));
86 
87  // Build importer from Domain Map to Column Map
88  col_importer =
92 
93  num_col_blocks = epetraCijk->numMyCols();
94  num_row_blocks = epetraCijk->numMyRows();
95  }
96 
97  input_block.resize(num_col_blocks);
98  result_block.resize(num_row_blocks);
99 }
100 
102 {
103 }
104 
105 void
109 {
110  block_ops = ops;
111  num_blocks = sg_basis->dimension() + 1;
112 }
113 
117 {
118  return block_ops;
119 }
120 
124 {
125  return block_ops;
126 }
127 
128 int
130 SetUseTranspose(bool UseTheTranspose)
131 {
132  useTranspose = UseTheTranspose;
133  for (int i=0; i<num_blocks; i++)
134  (*block_ops)[i].SetUseTranspose(useTranspose);
135 
136  return 0;
137 }
138 
139 int
141 Apply(const Epetra_MultiVector& Input, Epetra_MultiVector& Result) const
142 {
143  // We have to be careful if Input and Result are the same vector.
144  // If this is the case, the only possible solution is to make a copy
145  const Epetra_MultiVector *input = &Input;
146  bool made_copy = false;
147  if (Input.Values() == Result.Values() && !is_stoch_parallel) {
148  input = new Epetra_MultiVector(Input);
149  made_copy = true;
150  }
151 
152  // Initialize
153  Result.PutScalar(0.0);
154 
155  const Epetra_Map* input_base_map = domain_base_map.get();
156  const Epetra_Map* result_base_map = range_base_map.get();
157  if (useTranspose == true) {
158  input_base_map = range_base_map.get();
159  result_base_map = domain_base_map.get();
160  }
161 
162  // Allocate temporary storage
163  int m = Input.NumVectors();
164  if (useTranspose == false &&
165  (tmp == Teuchos::null || tmp->NumVectors() != m*max_num_mat_vec))
166  tmp = Teuchos::rcp(new Epetra_MultiVector(*result_base_map,
167  m*max_num_mat_vec));
168  else if (useTranspose == true &&
169  (tmp_trans == Teuchos::null ||
170  tmp_trans->NumVectors() != m*max_num_mat_vec))
171  tmp_trans = Teuchos::rcp(new Epetra_MultiVector(*result_base_map,
172  m*max_num_mat_vec));
173  Epetra_MultiVector *tmp_result;
174  if (useTranspose == false)
175  tmp_result = tmp.get();
176  else
177  tmp_result = tmp_trans.get();
178 
179  // Map input into column map
180  const Epetra_MultiVector *tmp_col;
181  if (!is_stoch_parallel)
182  tmp_col = input;
183  else {
184  if (useTranspose == false) {
185  if (input_col == Teuchos::null || input_col->NumVectors() != m)
186  input_col = Teuchos::rcp(new Epetra_MultiVector(*global_col_map, m));
187  input_col->Import(*input, *col_importer, Insert);
188  tmp_col = input_col.get();
189  }
190  else {
191  if (input_col_trans == Teuchos::null ||
192  input_col_trans->NumVectors() != m)
193  input_col_trans =
194  Teuchos::rcp(new Epetra_MultiVector(*global_col_map_trans, m));
195  input_col_trans->Import(*input, *col_importer_trans, Insert);
196  tmp_col = input_col_trans.get();
197  }
198  }
199 
200  // Extract blocks
201  EpetraExt::BlockMultiVector sg_input(View, *input_base_map, *tmp_col);
202  EpetraExt::BlockMultiVector sg_result(View, *result_base_map, Result);
203  for (int i=0; i<input_block.size(); i++)
204  input_block[i] = sg_input.GetBlock(i);
205  for (int i=0; i<result_block.size(); i++)
206  result_block[i] = sg_result.GetBlock(i);
207  int N = result_block[0]->MyLength();
208 
209  const Teuchos::Array<double>& norms = sg_basis->norm_squared();
210  int d = sg_basis->dimension();
211  Teuchos::Array<double> zero(d), one(d);
212  for(int j = 0; j<d; j++) {
213  zero[j] = 0.0;
214  one[j] = 1.0;
215  }
216  Teuchos::Array< double > phi_0(expansion_size), phi_1(expansion_size);
217  sg_basis->evaluateBases(zero, phi_0);
218  sg_basis->evaluateBases(one, phi_1);
219 
220  // k_begin and k_end are initialized in the constructor
221  for (Cijk_type::k_iterator k_it=k_begin; k_it!=k_end; ++k_it) {
222  Cijk_type::kj_iterator j_begin = Cijk->j_begin(k_it);
223  Cijk_type::kj_iterator j_end = Cijk->j_end(k_it);
224  int k = index(k_it);
225  int nj = Cijk->num_j(k_it);
226  if (nj > 0) {
227  Teuchos::Array<double*> j_ptr(nj*m);
228  Teuchos::Array<int> mj_indices(nj*m);
229  int l = 0;
230  for (Cijk_type::kj_iterator j_it = j_begin; j_it != j_end; ++j_it) {
231  int j = index(j_it);
232  for (int mm=0; mm<m; mm++) {
233  j_ptr[l*m+mm] = input_block[j]->Values()+mm*N;
234  mj_indices[l*m+mm] = l*m+mm;
235  }
236  l++;
237  }
238  Epetra_MultiVector input_tmp(View, *input_base_map, &j_ptr[0], nj*m);
239  Epetra_MultiVector result_tmp(View, *tmp_result, &mj_indices[0], nj*m);
240  (*block_ops)[k].Apply(input_tmp, result_tmp);
241  l = 0;
242  for (Cijk_type::kj_iterator j_it = j_begin; j_it != j_end; ++j_it) {
243  int j = index(j_it);
244  int j_gid = epetraCijk->GCID(j);
245  for (Cijk_type::kji_iterator i_it = Cijk->i_begin(j_it);
246  i_it != Cijk->i_end(j_it); ++i_it) {
247  int i = index(i_it);
248  int i_gid = epetraCijk->GRID(i);
249  double c = value(i_it);
250  if (k == 0)
251  c /= phi_0[0];
252  else {
253  c /= phi_1[k];
254  if (i_gid == j_gid)
255  c -= phi_0[k]/(phi_1[k]*phi_0[0])*norms[i_gid];
256  }
257  if (scale_op) {
258  if (useTranspose)
259  c /= norms[j_gid];
260  else
261  c /= norms[i_gid];
262  }
263  for (int mm=0; mm<m; mm++)
264  (*result_block[i])(mm)->Update(c, *result_tmp(l*m+mm), 1.0);
265  }
266  l++;
267  }
268  }
269  }
270 
271  // Destroy blocks
272  for (int i=0; i<input_block.size(); i++)
273  input_block[i] = Teuchos::null;
274  for (int i=0; i<result_block.size(); i++)
275  result_block[i] = Teuchos::null;
276 
277  if (made_copy)
278  delete input;
279 
280  return 0;
281 }
282 
283 int
286 {
287  throw "KLMatrixFreeOperator::ApplyInverse not defined!";
288  return -1;
289 }
290 
291 double
293 NormInf() const
294 {
295  return 1.0;
296 }
297 
298 
299 const char*
301 Label () const
302 {
303  return const_cast<char*>(label.c_str());
304 }
305 
306 bool
309 {
310  return useTranspose;
311 }
312 
313 bool
315 HasNormInf() const
316 {
317  return false;
318 }
319 
320 const Epetra_Comm &
322 Comm() const
323 {
324  return *sg_comm;
325 }
326 const Epetra_Map&
329 {
330  if (useTranspose)
331  return *range_sg_map;
332  return *domain_sg_map;
333 }
334 
335 const Epetra_Map&
338 {
339  if (useTranspose)
340  return *domain_sg_map;
341  return *range_sg_map;
342 }
Teuchos::RCP< const Epetra_Map > range_sg_map
Stores range SG map.
virtual int Apply(const Epetra_MultiVector &Input, Epetra_MultiVector &Result) const
Returns the result of a Epetra_Operator applied to a Epetra_MultiVector Input in Result as described ...
Cijk_type::k_iterator k_end
Ending k iterator.
Teuchos::RCP< Epetra_Map > global_col_map_trans
Stores operator column SG map for transpose.
bool scale_op
Flag indicating whether operator be scaled with &lt;^2&gt;
virtual double NormInf() const
Returns an approximate infinity norm of the operator matrix.
Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > epetraCijk
Stores Epetra Cijk tensor.
bool is_stoch_parallel
Whether we have parallelism over stochastic blocks.
T & get(ParameterList &l, const std::string &name)
Teuchos::RCP< const EpetraExt::MultiComm > sg_comm
Stores SG parallel communicator.
virtual int SetUseTranspose(bool UseTranspose)
Set to true if the transpose of the operator is requested.
virtual Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > getSGPolynomial()
Get SG polynomial.
KLMatrixFreeOperator(const Teuchos::RCP< const EpetraExt::MultiComm > &sg_comm, const Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > &sg_basis, const Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > &epetraCijk, const Teuchos::RCP< const Epetra_Map > &domain_base_map, const Teuchos::RCP< const Epetra_Map > &range_base_map, const Teuchos::RCP< const Epetra_Map > &domain_sg_map, const Teuchos::RCP< const Epetra_Map > &range_sg_map, const Teuchos::RCP< Teuchos::ParameterList > &params)
Constructor.
virtual ordinal_type dimension() const =0
Return dimension of basis.
virtual const char * Label() const
Returns a character std::string describing the operator.
int expansion_size
Number of terms in expansion.
Teuchos::RCP< const Epetra_Map > range_base_map
Stores range base map.
Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > sg_basis
Stochastic Galerking basis.
Bi-directional iterator for traversing a sparse array.
Teuchos::Array< Teuchos::RCP< const Epetra_MultiVector > > input_block
MultiVectors for each block for Apply() input.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< const Epetra_Map > domain_sg_map
Stores domain SG map.
virtual bool UseTranspose() const
Returns the current UseTranspose setting.
Cijk_type::k_iterator k_begin
Starting k iterator.
int max_num_mat_vec
Maximum number of matvecs in Apply.
Teuchos::RCP< const Cijk_type > Cijk
Stores triple product tensor.
Teuchos::RCP< const Epetra_BlockMap > stoch_col_map
Stores stochastic part of column map.
Teuchos::RCP< Epetra_Import > col_importer
Importer from domain map to column map.
virtual const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this matrix operator. ...
virtual bool HasNormInf() const
Returns true if the this object can provide an approximate Inf-norm, false otherwise.
bool include_mean
Flag indicating whether to include mean term.
virtual const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this matrix operator.
Teuchos::RCP< Epetra_Import > col_importer_trans
Importer from range map to column map.
Teuchos::RCP< Epetra_Map > global_col_map
Stores operator column SG map.
Teuchos::RCP< const Epetra_Map > domain_base_map
Stores domain base map.
virtual void setupOperator(const Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > &poly)
Setup operator.
virtual int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of the inverse of the operator applied to a Epetra_MultiVector Input in Result as ...
Teuchos::Array< Teuchos::RCP< Epetra_MultiVector > > result_block
MultiVectors for each block for Apply() result.
virtual const Epetra_Comm & Comm() const
Returns a reference to the Epetra_Comm communicator associated with this operator.