trilinos_epetra_vector.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 // C++ includes
19 #include <limits>
20 
21 // Local Includes
23 
24 #ifdef LIBMESH_TRILINOS_HAVE_EPETRA
25 
27 #include "libmesh/dense_vector.h"
28 #include "libmesh/parallel.h"
30 #include "libmesh/utility.h"
31 
32 // Trilinos Includes
34 #include <Epetra_LocalMap.h>
35 #include <Epetra_Comm.h>
36 #include <Epetra_Map.h>
37 #include <Epetra_BlockMap.h>
38 #include <Epetra_Import.h>
39 #include <Epetra_Export.h>
40 #include <Epetra_Util.h>
41 #include <Epetra_IntSerialDenseVector.h>
42 #include <Epetra_SerialDenseVector.h>
43 #include <Epetra_Vector.h>
45 
46 namespace libMesh
47 {
48 
49 template <typename T>
51 {
52  libmesh_assert(this->closed());
53 
54  const unsigned int nl = _vec->MyLength();
55 
56  T sum=0.0;
57 
58  T * values = _vec->Values();
59 
60  for (unsigned int i=0; i<nl; i++)
61  sum += values[i];
62 
63  this->comm().sum(sum);
64 
65  return sum;
66 }
67 
68 template <typename T>
70 {
71  libmesh_assert(this->closed());
72 
73  Real value;
74 
75  _vec->Norm1(&value);
76 
77  return value;
78 }
79 
80 template <typename T>
82 {
83  libmesh_assert(this->closed());
84 
85  Real value;
86 
87  _vec->Norm2(&value);
88 
89  return value;
90 }
91 
92 template <typename T>
94 {
95  libmesh_assert(this->closed());
96 
97  Real value;
98 
99  _vec->NormInf(&value);
100 
101  return value;
102 }
103 
104 template <typename T>
107 {
108  libmesh_assert(this->closed());
109 
110  this->add(1., v);
111 
112  return *this;
113 }
114 
115 
116 
117 template <typename T>
120 {
121  libmesh_assert(this->closed());
122 
123  this->add(-1., v);
124 
125  return *this;
126 }
127 
128 
129 template <typename T>
132 {
133  libmesh_assert(this->closed());
134  libmesh_assert_equal_to(size(), v.size());
135 
136  const EpetraVector<T> & v_vec = cast_ref<const EpetraVector<T> &>(v);
137 
138  _vec->ReciprocalMultiply(1.0, *v_vec._vec, *_vec, 0.0);
139 
140  return *this;
141 }
142 
143 
144 
145 
146 template <typename T>
147 void EpetraVector<T>::set (const numeric_index_type i_in, const T value_in)
148 {
149  int i = static_cast<int> (i_in);
150  T value = value_in;
151 
152  libmesh_assert_less (i_in, this->size());
153 
154  ReplaceGlobalValues(1, &i, &value);
155 
156  this->_is_closed = false;
157 }
158 
159 
160 
161 template <typename T>
163 {
164  // The Epetra::reciprocal() function takes a constant reference to *another* vector,
165  // and fills _vec with its reciprocal. Does that mean we can't pass *_vec as the
166  // argument?
167  // _vec->reciprocal( *_vec );
168 
169  // Alternatively, compute the reciprocal by hand... see also the add(T) member that does this...
170  const unsigned int nl = _vec->MyLength();
171 
172  T * values = _vec->Values();
173 
174  for (unsigned int i=0; i<nl; i++)
175  {
176  // Don't divide by zero (maybe only check this in debug mode?)
177  if (std::abs(values[i]) < std::numeric_limits<T>::min())
178  libmesh_error_msg("Error, divide by zero in DistributedVector<T>::reciprocal()!");
179 
180  values[i] = 1. / values[i];
181  }
182 
183  // Leave the vector in a closed state...
184  this->close();
185 }
186 
187 
188 
189 template <typename T>
191 {
192  // EPetra is real, rendering this a no-op.
193 }
194 
195 
196 
197 template <typename T>
198 void EpetraVector<T>::add (const numeric_index_type i_in, const T value_in)
199 {
200  int i = static_cast<int> (i_in);
201  T value = value_in;
202 
203  libmesh_assert_less (i_in, this->size());
204 
205  SumIntoGlobalValues(1, &i, &value);
206 
207  this->_is_closed = false;
208 }
209 
210 
211 
212 template <typename T>
213 void EpetraVector<T>::add_vector (const T * v,
214  const std::vector<numeric_index_type> & dof_indices)
215 {
216  libmesh_assert_equal_to (sizeof(numeric_index_type), sizeof(int));
217 
218  SumIntoGlobalValues (cast_int<numeric_index_type>(dof_indices.size()),
219  numeric_trilinos_cast(dof_indices.data()),
220  const_cast<T *>(v));
221 }
222 
223 
224 
225 // TODO: fill this in after creating an EpetraMatrix
226 template <typename T>
228  const SparseMatrix<T> & A_in)
229 {
230  const EpetraVector<T> * v = cast_ptr<const EpetraVector<T> *>(&v_in);
231  const EpetraMatrix<T> * A = cast_ptr<const EpetraMatrix<T> *>(&A_in);
232 
233  // FIXME - does Trilinos let us do this *without* memory allocation?
234  std::unique_ptr<NumericVector<T>> temp = v->zero_clone();
235  EpetraVector<T> * temp_v = cast_ptr<EpetraVector<T> *>(temp.get());
236  A->mat()->Multiply(false, *v->_vec, *temp_v->_vec);
237  *this += *temp;
238 }
239 
240 
241 
242 // TODO: fill this in after creating an EpetraMatrix
243 template <typename T>
245  const SparseMatrix<T> & /* A_in */)
246 {
247  libmesh_not_implemented();
248 }
249 
250 
251 
252 template <typename T>
253 void EpetraVector<T>::add (const T v_in)
254 {
255  const unsigned int nl = _vec->MyLength();
256 
257  T * values = _vec->Values();
258 
259  for (unsigned int i=0; i<nl; i++)
260  values[i]+=v_in;
261 
262  this->_is_closed = false;
263 }
264 
265 
266 template <typename T>
268 {
269  this->add (1., v);
270 }
271 
272 
273 template <typename T>
274 void EpetraVector<T>::add (const T a_in, const NumericVector<T> & v_in)
275 {
276  const EpetraVector<T> * v = cast_ptr<const EpetraVector<T> *>(&v_in);
277 
278  libmesh_assert_equal_to (this->size(), v->size());
279 
280  _vec->Update(a_in,*v->_vec, 1.);
281 }
282 
283 
284 
285 template <typename T>
286 void EpetraVector<T>::insert (const T * v,
287  const std::vector<numeric_index_type> & dof_indices)
288 {
289  libmesh_assert_equal_to (sizeof(numeric_index_type), sizeof(int));
290 
291  ReplaceGlobalValues (cast_int<numeric_index_type>(dof_indices.size()),
292  numeric_trilinos_cast(dof_indices.data()),
293  const_cast<T *>(v));
294 }
295 
296 
297 
298 template <typename T>
299 void EpetraVector<T>::scale (const T factor_in)
300 {
301  _vec->Scale(factor_in);
302 }
303 
304 template <typename T>
306 {
307  _vec->Abs(*_vec);
308 }
309 
310 
311 template <typename T>
313 {
314  const EpetraVector<T> * v = cast_ptr<const EpetraVector<T> *>(&v_in);
315 
316  T result=0.0;
317 
318  _vec->Dot(*v->_vec, &result);
319 
320  return result;
321 }
322 
323 
324 template <typename T>
326  const NumericVector<T> & vec2)
327 {
328  const EpetraVector<T> * v1 = cast_ptr<const EpetraVector<T> *>(&vec1);
329  const EpetraVector<T> * v2 = cast_ptr<const EpetraVector<T> *>(&vec2);
330 
331  _vec->Multiply(1.0, *v1->_vec, *v2->_vec, 0.0);
332 }
333 
334 
335 template <typename T>
338 {
339  _vec->PutScalar(s_in);
340 
341  return *this;
342 }
343 
344 
345 
346 template <typename T>
349 {
350  // This function could be implemented in terms of the copy
351  // assignment operator (see other NumericVector subclasses) but that
352  // function is currently deleted so calling this function is an error.
353  // const EpetraVector<T> * v = cast_ptr<const EpetraVector<T> *>(&v_in);
354  // *this = *v;
355  libmesh_not_implemented();
356  return *this;
357 }
358 
359 
360 
361 template <typename T>
363 EpetraVector<T>::operator = (const std::vector<T> & v)
364 {
365  T * values = _vec->Values();
366 
371  if (this->size() == v.size())
372  {
373  const unsigned int nl=this->local_size();
374  const unsigned int fli=this->first_local_index();
375 
376  for (unsigned int i=0;i<nl;i++)
377  values[i]=v[fli+i];
378  }
379 
384  else
385  {
386  libmesh_assert_equal_to (v.size(), this->local_size());
387 
388  const unsigned int nl=this->local_size();
389 
390  for (unsigned int i=0;i<nl;i++)
391  values[i]=v[i];
392  }
393 
394  return *this;
395 }
396 
397 
398 
399 template <typename T>
401 {
402  EpetraVector<T> * v_local = cast_ptr<EpetraVector<T> *>(&v_local_in);
403 
404  Epetra_Map rootMap = Epetra_Util::Create_Root_Map( *_map, -1);
405  v_local->_vec->ReplaceMap(rootMap);
406 
407  Epetra_Import importer(v_local->_vec->Map(), *_map);
408  v_local->_vec->Import(*_vec, importer, Insert);
409 }
410 
411 
412 
413 template <typename T>
415  const std::vector<numeric_index_type> & /* send_list */) const
416 {
417  // TODO: optimize to sync only the send list values
418  this->localize(v_local_in);
419 
420  // EpetraVector<T> * v_local =
421  // cast_ptr<EpetraVector<T> *>(&v_local_in);
422 
423  // libmesh_assert(this->_map.get());
424  // libmesh_assert(v_local->_map.get());
425  // libmesh_assert_equal_to (v_local->local_size(), this->size());
426  // libmesh_assert_less_equal (send_list.size(), v_local->size());
427 
428  // Epetra_Import importer (*v_local->_map, *this->_map);
429 
430  // v_local->_vec->Import (*this->_vec, importer, Insert);
431 }
432 
433 
434 
435 template <typename T>
436 void EpetraVector<T>::localize (std::vector<T> & v_local,
437  const std::vector<numeric_index_type> & indices) const
438 {
439  // Create a "replicated" map for importing values. This is
440  // equivalent to creating a general Epetra_Map with
441  // NumGlobalElements == NumMyElements.
442  Epetra_LocalMap import_map(static_cast<int>(indices.size()),
443  /*IndexBase=*/0,
444  _map->Comm());
445 
446  // Get a pointer to the list of global elements for the map, and set
447  // all the values from indices.
448  int * import_map_global_elements = import_map.MyGlobalElements();
449  for (int i=0; i<import_map.NumMyElements(); ++i)
450  import_map_global_elements[i] = indices[i];
451 
452  // Create a new EpetraVector to import values into.
453  Epetra_Vector import_vector(import_map);
454 
455  // Set up an "Import" object which associates the two maps.
456  Epetra_Import import_object(import_map, *_map);
457 
458  // Import the values
459  import_vector.Import(*_vec, import_object, Insert);
460 
461  // Get a pointer to the imported values array and the length of the
462  // array.
463  T * values = import_vector.Values();
464  int import_vector_length = import_vector.MyLength();
465 
466  // Copy the imported values into v_local
467  v_local.resize(import_vector_length);
468  for (int i=0; i<import_vector_length; ++i)
469  v_local[i] = values[i];
470 }
471 
472 
473 
474 template <typename T>
475 void EpetraVector<T>::localize (const numeric_index_type first_local_idx,
476  const numeric_index_type last_local_idx,
477  const std::vector<numeric_index_type> & send_list)
478 {
479  // Only good for serial vectors.
480  libmesh_assert_equal_to (this->size(), this->local_size());
481  libmesh_assert_greater (last_local_idx, first_local_idx);
482  libmesh_assert_less_equal (send_list.size(), this->size());
483  libmesh_assert_less (last_local_idx, this->size());
484 
485  const unsigned int my_size = this->size();
486  const unsigned int my_local_size = (last_local_idx - first_local_idx + 1);
487 
488  // Don't bother for serial cases
489  if ((first_local_idx == 0) &&
490  (my_local_size == my_size))
491  return;
492 
493  // Build a parallel vector, initialize it with the local
494  // parts of (*this)
495  EpetraVector<T> parallel_vec(this->comm(), PARALLEL);
496 
497  parallel_vec.init (my_size, my_local_size, true, PARALLEL);
498 
499  // Copy part of *this into the parallel_vec
500  for (numeric_index_type i=first_local_idx; i<=last_local_idx; i++)
501  parallel_vec.set(i,this->el(i));
502 
503  // localize like normal
504  parallel_vec.close();
505  parallel_vec.localize (*this, send_list);
506  this->close();
507 }
508 
509 
510 
511 template <typename T>
512 void EpetraVector<T>::localize (std::vector<T> & v_local) const
513 {
514  // This function must be run on all processors at once
515  parallel_object_only();
516 
517  const unsigned int n = this->size();
518  const unsigned int nl = this->local_size();
519 
520  libmesh_assert(this->_vec);
521 
522  v_local.clear();
523  v_local.reserve(n);
524 
525  // build up my local part
526  for (unsigned int i=0; i<nl; i++)
527  v_local.push_back((*this->_vec)[i]);
528 
529  this->comm().allgather (v_local);
530 }
531 
532 
533 
534 template <typename T>
535 void EpetraVector<T>::localize_to_one (std::vector<T> & v_local,
536  const processor_id_type pid) const
537 {
538  // This function must be run on all processors at once
539  parallel_object_only();
540 
541  const unsigned int n = this->size();
542  const unsigned int nl = this->local_size();
543 
544  libmesh_assert_less (pid, this->n_processors());
545  libmesh_assert(this->_vec);
546 
547  v_local.clear();
548  v_local.reserve(n);
549 
550 
551  // build up my local part
552  for (unsigned int i=0; i<nl; i++)
553  v_local.push_back((*this->_vec)[i]);
554 
555  this->comm().gather (pid, v_local);
556 }
557 
558 
559 
560 template <typename T>
562  const std::vector<numeric_index_type> & /* rows */) const
563 {
564  libmesh_not_implemented();
565 }
566 
567 
568 /*********************************************************************
569  * The following were copied (and slightly modified) from
570  * Epetra_FEVector.h in order to allow us to use a standard
571  * Epetra_Vector... which is more compatible with other Trilinos
572  * packages such as NOX. All of this code is originally under LGPL
573  *********************************************************************/
574 
575 //----------------------------------------------------------------------------
576 template <typename T>
578  const int * GIDs,
579  const double * values)
580 {
581  return( inputValues( numIDs, GIDs, values, true) );
582 }
583 
584 //----------------------------------------------------------------------------
585 template <typename T>
586 int EpetraVector<T>::SumIntoGlobalValues(const Epetra_IntSerialDenseVector & GIDs,
587  const Epetra_SerialDenseVector & values)
588 {
589  if (GIDs.Length() != values.Length()) {
590  return(-1);
591  }
592 
593  return( inputValues( GIDs.Length(), GIDs.Values(), values.Values(), true) );
594 }
595 
596 //----------------------------------------------------------------------------
597 template <typename T>
599  const int * GIDs,
600  const int * numValuesPerID,
601  const double * values)
602 {
603  return( inputValues( numIDs, GIDs, numValuesPerID, values, true) );
604 }
605 
606 //----------------------------------------------------------------------------
607 template <typename T>
609  const int * GIDs,
610  const double * values)
611 {
612  return( inputValues( numIDs, GIDs, values, false) );
613 }
614 
615 //----------------------------------------------------------------------------
616 template <typename T>
617 int EpetraVector<T>::ReplaceGlobalValues(const Epetra_IntSerialDenseVector & GIDs,
618  const Epetra_SerialDenseVector & values)
619 {
620  if (GIDs.Length() != values.Length()) {
621  return(-1);
622  }
623 
624  return( inputValues( GIDs.Length(), GIDs.Values(), values.Values(), false) );
625 }
626 
627 //----------------------------------------------------------------------------
628 template <typename T>
630  const int * GIDs,
631  const int * numValuesPerID,
632  const double * values)
633 {
634  return( inputValues( numIDs, GIDs, numValuesPerID, values, false) );
635 }
636 
637 //----------------------------------------------------------------------------
638 template <typename T>
640  const int * GIDs,
641  const double * values,
642  bool accumulate)
643 {
644  if (accumulate) {
645  libmesh_assert(last_edit == 0 || last_edit == 2);
646  last_edit = 2;
647  } else {
648  libmesh_assert(last_edit == 0 || last_edit == 1);
649  last_edit = 1;
650  }
651 
652  //Important note!! This method assumes that there is only 1 point
653  //associated with each element.
654 
655  for (int i=0; i<numIDs; ++i) {
656  if (_vec->Map().MyGID(GIDs[i])) {
657  if (accumulate) {
658  _vec->SumIntoGlobalValue(GIDs[i], 0, 0, values[i]);
659  }
660  else {
661  _vec->ReplaceGlobalValue(GIDs[i], 0, 0, values[i]);
662  }
663  }
664  else {
665  if (!ignoreNonLocalEntries_) {
666  EPETRA_CHK_ERR( inputNonlocalValue(GIDs[i], values[i], accumulate) );
667  }
668  }
669  }
670 
671  return(0);
672 }
673 
674 //----------------------------------------------------------------------------
675 template <typename T>
677  const int * GIDs,
678  const int * numValuesPerID,
679  const double * values,
680  bool accumulate)
681 {
682  if (accumulate) {
683  libmesh_assert(last_edit == 0 || last_edit == 2);
684  last_edit = 2;
685  } else {
686  libmesh_assert(last_edit == 0 || last_edit == 1);
687  last_edit = 1;
688  }
689 
690  int offset=0;
691  for (int i=0; i<numIDs; ++i) {
692  int numValues = numValuesPerID[i];
693  if (_vec->Map().MyGID(GIDs[i])) {
694  if (accumulate) {
695  for (int j=0; j<numValues; ++j) {
696  _vec->SumIntoGlobalValue(GIDs[i], j, 0, values[offset+j]);
697  }
698  }
699  else {
700  for (int j=0; j<numValues; ++j) {
701  _vec->ReplaceGlobalValue(GIDs[i], j, 0, values[offset+j]);
702  }
703  }
704  }
705  else {
706  if (!ignoreNonLocalEntries_) {
707  EPETRA_CHK_ERR( inputNonlocalValues(GIDs[i], numValues,
708  &(values[offset]), accumulate) );
709  }
710  }
711  offset += numValues;
712  }
713 
714  return(0);
715 }
716 
717 //----------------------------------------------------------------------------
718 template <typename T>
719 int EpetraVector<T>::inputNonlocalValue(int GID, double value, bool accumulate)
720 {
721  int insertPoint = -1;
722 
723  //find offset of GID in nonlocalIDs_
724  int offset = Epetra_Util_binary_search(GID, nonlocalIDs_, numNonlocalIDs_,
725  insertPoint);
726  if (offset >= 0) {
727  //if offset >= 0
728  // put value in nonlocalCoefs_[offset][0]
729 
730  if (accumulate) {
731  nonlocalCoefs_[offset][0] += value;
732  }
733  else {
734  nonlocalCoefs_[offset][0] = value;
735  }
736  }
737  else {
738  //else
739  // insert GID in nonlocalIDs_
740  // insert 1 in nonlocalElementSize_
741  // insert value in nonlocalCoefs_
742 
743  int tmp1 = numNonlocalIDs_;
744  int tmp2 = allocatedNonlocalLength_;
745  int tmp3 = allocatedNonlocalLength_;
746  EPETRA_CHK_ERR( Epetra_Util_insert(GID, insertPoint, nonlocalIDs_,
747  tmp1, tmp2) );
748  --tmp1;
749  EPETRA_CHK_ERR( Epetra_Util_insert(1, insertPoint, nonlocalElementSize_,
750  tmp1, tmp3) );
751  double * values = new double[1];
752  values[0] = value;
753  EPETRA_CHK_ERR( Epetra_Util_insert(values, insertPoint, nonlocalCoefs_,
754  numNonlocalIDs_, allocatedNonlocalLength_) );
755  }
756 
757  return(0);
758 }
759 
760 //----------------------------------------------------------------------------
761 template <typename T>
763  int numValues,
764  const double * values,
765  bool accumulate)
766 {
767  int insertPoint = -1;
768 
769  //find offset of GID in nonlocalIDs_
770  int offset = Epetra_Util_binary_search(GID, nonlocalIDs_, numNonlocalIDs_,
771  insertPoint);
772  if (offset >= 0) {
773  //if offset >= 0
774  // put value in nonlocalCoefs_[offset][0]
775 
776  if (numValues != nonlocalElementSize_[offset]) {
777  libMesh::err << "Epetra_FEVector ERROR: block-size for GID " << GID << " is "
778  << numValues<<" which doesn't match previously set block-size of "
779  << nonlocalElementSize_[offset] << std::endl;
780  return(-1);
781  }
782 
783  if (accumulate) {
784  for (int j=0; j<numValues; ++j) {
785  nonlocalCoefs_[offset][j] += values[j];
786  }
787  }
788  else {
789  for (int j=0; j<numValues; ++j) {
790  nonlocalCoefs_[offset][j] = values[j];
791  }
792  }
793  }
794  else {
795  //else
796  // insert GID in nonlocalIDs_
797  // insert numValues in nonlocalElementSize_
798  // insert values in nonlocalCoefs_
799 
800  int tmp1 = numNonlocalIDs_;
801  int tmp2 = allocatedNonlocalLength_;
802  int tmp3 = allocatedNonlocalLength_;
803  EPETRA_CHK_ERR( Epetra_Util_insert(GID, insertPoint, nonlocalIDs_,
804  tmp1, tmp2) );
805  --tmp1;
806  EPETRA_CHK_ERR( Epetra_Util_insert(numValues, insertPoint, nonlocalElementSize_,
807  tmp1, tmp3) );
808  double * newvalues = new double[numValues];
809  for (int j=0; j<numValues; ++j) {
810  newvalues[j] = values[j];
811  }
812  EPETRA_CHK_ERR( Epetra_Util_insert(newvalues, insertPoint, nonlocalCoefs_,
813  numNonlocalIDs_, allocatedNonlocalLength_) );
814  }
815 
816  return(0);
817 }
818 
819 //----------------------------------------------------------------------------
820 template <typename T>
821 int EpetraVector<T>::GlobalAssemble(Epetra_CombineMode mode)
822 {
823  //In this method we need to gather all the non-local (overlapping) data
824  //that's been input on each processor, into the (probably) non-overlapping
825  //distribution defined by the map that 'this' vector was constructed with.
826 
827  //We don't need to do anything if there's only one processor or if
828  //ignoreNonLocalEntries_ is true.
829  if (_vec->Map().Comm().NumProc() < 2 || ignoreNonLocalEntries_) {
830  return(0);
831  }
832 
833 
834 
835  //First build a map that describes the data in nonlocalIDs_/nonlocalCoefs_.
836  //We'll use the arbitrary distribution constructor of Map.
837 
838  Epetra_BlockMap sourceMap(-1, numNonlocalIDs_,
839  nonlocalIDs_, nonlocalElementSize_,
840  _vec->Map().IndexBase(), _vec->Map().Comm());
841 
842  //Now build a vector to hold our nonlocalCoefs_, and to act as the source-
843  //vector for our import operation.
844  Epetra_MultiVector nonlocalVector(sourceMap, 1);
845 
846  int i,j;
847  for (i=0; i<numNonlocalIDs_; ++i) {
848  for (j=0; j<nonlocalElementSize_[i]; ++j) {
849  nonlocalVector.ReplaceGlobalValue(nonlocalIDs_[i], j, 0,
850  nonlocalCoefs_[i][j]);
851  }
852  }
853 
854  Epetra_Export exporter(sourceMap, _vec->Map());
855 
856  EPETRA_CHK_ERR( _vec->Export(nonlocalVector, exporter, mode) );
857 
858  destroyNonlocalData();
859 
860  return(0);
861 }
862 
863 
864 //----------------------------------------------------------------------------
865 template <typename T>
867 {
868  (*_vec) = *(source._vec);
869 
870  destroyNonlocalData();
871 
872  if (source.allocatedNonlocalLength_ > 0) {
873  allocatedNonlocalLength_ = source.allocatedNonlocalLength_;
874  numNonlocalIDs_ = source.numNonlocalIDs_;
875  nonlocalIDs_ = new int[allocatedNonlocalLength_];
876  nonlocalElementSize_ = new int[allocatedNonlocalLength_];
877  nonlocalCoefs_ = new double *[allocatedNonlocalLength_];
878  for (int i=0; i<numNonlocalIDs_; ++i) {
879  int elemSize = source.nonlocalElementSize_[i];
880  nonlocalCoefs_[i] = new double[elemSize];
881  nonlocalIDs_[i] = source.nonlocalIDs_[i];
882  nonlocalElementSize_[i] = elemSize;
883  for (int j=0; j<elemSize; ++j) {
884  nonlocalCoefs_[i][j] = source.nonlocalCoefs_[i][j];
885  }
886  }
887  }
888 }
889 
890 
891 //----------------------------------------------------------------------------
892 template <typename T>
894 {
895  if (allocatedNonlocalLength_ > 0) {
896  delete [] nonlocalIDs_;
897  delete [] nonlocalElementSize_;
898  nonlocalIDs_ = nullptr;
899  nonlocalElementSize_ = nullptr;
900  for (int i=0; i<numNonlocalIDs_; ++i) {
901  delete [] nonlocalCoefs_[i];
902  }
903  delete [] nonlocalCoefs_;
904  nonlocalCoefs_ = nullptr;
905  numNonlocalIDs_ = 0;
906  allocatedNonlocalLength_ = 0;
907  }
908  return;
909 }
910 
911 
912 //------------------------------------------------------------------
913 // Explicit instantiations
914 template class EpetraVector<Number>;
915 
916 } // namespace libMesh
917 
918 #endif // LIBMESH_TRILINOS_HAVE_EPETRA
virtual NumericVector< T > & operator/=(const NumericVector< T > &v) override
bool closed()
Definition: libmesh.C:265
double abs(double a)
virtual void pointwise_mult(const NumericVector< T > &vec1, const NumericVector< T > &vec2) override
virtual void reciprocal() override
virtual void add(const numeric_index_type i, const T value) override
virtual void close() override
virtual std::unique_ptr< NumericVector< T > > zero_clone() const override
virtual numeric_index_type size() const =0
virtual void localize_to_one(std::vector< T > &v_local, const processor_id_type proc_id=0) const override
virtual void set(const numeric_index_type i, const T value) override
uint8_t processor_id_type
Definition: id_types.h:99
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
Definition: diff_context.h:40
virtual void insert(const T *v, const std::vector< numeric_index_type > &dof_indices) override
int inputNonlocalValues(int GID, int numValues, const double *values, bool accumulate)
virtual void conjugate() override
virtual T sum() const override
int GlobalAssemble(Epetra_CombineMode mode=Add)
virtual numeric_index_type size() const override
virtual void add_vector_transpose(const NumericVector< T > &v, const SparseMatrix< T > &A) override
int ReplaceGlobalValues(int numIDs, const int *GIDs, const double *values)
int inputValues(int numIDs, const int *GIDs, const double *values, bool accumulate)
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual NumericVector< T > & operator+=(const NumericVector< T > &v) override
virtual Real l1_norm() const override
virtual void abs() override
virtual NumericVector< T > & operator-=(const NumericVector< T > &v) override
int inputNonlocalValue(int GID, double value, bool accumulate)
virtual Real linfty_norm() const override
OStreamProxy err(std::cerr)
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType type=AUTOMATIC) override
virtual void localize(std::vector< T > &v_local) const override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual T dot(const NumericVector< T > &v) const override
int SumIntoGlobalValues(int numIDs, const int *GIDs, const double *values)
static PetscErrorCode Mat * A
virtual Real l2_norm() const override
static const bool value
Definition: xdr_io.C:109
virtual void create_subvector(NumericVector< T > &subvector, const std::vector< numeric_index_type > &rows) const override
void FEoperatorequals(const EpetraVector &source)
int * numeric_trilinos_cast(const numeric_index_type *p)
virtual void scale(const T factor) override
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices) override
long double min(long double a, double b)
EpetraVector & operator=(const EpetraVector &)=delete