Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Epetra_FECrsMatrix.cpp
Go to the documentation of this file.
1 
2 //@HEADER
3 // ************************************************************************
4 //
5 // Epetra: Linear Algebra Services Package
6 // Copyright 2011 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
39 //
40 // ************************************************************************
41 //@HEADER
42 
43 #include <Epetra_ConfigDefs.h>
44 #include <Epetra_FECrsMatrix.h>
47 #include <Epetra_FECrsGraph.h>
48 #include <Epetra_Export.h>
49 #include <Epetra_Comm.h>
50 #include <Epetra_Map.h>
51 #include <Epetra_Util.h>
52 
53 //----------------------------------------------------------------------------
55  const Epetra_Map& rowMap,
56  int* NumEntriesPerRow,
57  bool ignoreNonLocalEntries)
58  : Epetra_CrsMatrix(CV, rowMap, NumEntriesPerRow),
59  myFirstRow_(0),
60  myNumRows_(0),
61  ignoreNonLocalEntries_(ignoreNonLocalEntries),
62 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
65 #endif
66 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
69 #endif
71  workData_(128),
72  useNonlocalMatrix_ (false),
73  nonlocalMatrix_ (NULL),
74  sourceMap_(NULL),
75  colMap_(NULL),
76  exporter_(NULL),
77  tempMat_(NULL)
78 {
79  myFirstRow_ = rowMap.MinMyGID64();
80  myNumRows_ = rowMap.NumMyElements();
81 }
82 
83 //----------------------------------------------------------------------------
85  const Epetra_Map& rowMap,
86  int NumEntriesPerRow,
87  bool ignoreNonLocalEntries)
88  : Epetra_CrsMatrix(CV, rowMap, NumEntriesPerRow),
89  myFirstRow_(0),
90  myNumRows_(0),
91  ignoreNonLocalEntries_(ignoreNonLocalEntries),
92 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
95 #endif
96 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
99 #endif
100  nonlocalCoefs_(),
101  workData_(128),
102  useNonlocalMatrix_ (false),
103  nonlocalMatrix_ (NULL),
104  sourceMap_(NULL),
105  colMap_(NULL),
106  exporter_(NULL),
107  tempMat_(NULL)
108 {
109  myFirstRow_ = rowMap.MinMyGID64();
110  myNumRows_ = rowMap.NumMyElements();
111 }
112 
113 //----------------------------------------------------------------------------
115  const Epetra_Map& rowMap,
116  const Epetra_Map& colMap,
117  int* NumEntriesPerRow,
118  bool ignoreNonLocalEntries)
119  : Epetra_CrsMatrix(CV, rowMap, colMap, NumEntriesPerRow),
120  myFirstRow_(0),
121  myNumRows_(0),
122  ignoreNonLocalEntries_(ignoreNonLocalEntries),
123 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
126 #endif
127 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
130 #endif
131  nonlocalCoefs_(),
132  workData_(128),
133  useNonlocalMatrix_ (false),
134  nonlocalMatrix_ (NULL),
135  sourceMap_(NULL),
136  colMap_(NULL),
137  exporter_(NULL),
138  tempMat_(NULL)
139 {
140  myFirstRow_ = rowMap.MinMyGID64();
141  myNumRows_ = rowMap.NumMyElements();
142 }
143 
144 //----------------------------------------------------------------------------
146  const Epetra_Map& rowMap,
147  const Epetra_Map& colMap,
148  int NumEntriesPerRow,
149  bool ignoreNonLocalEntries)
150  : Epetra_CrsMatrix(CV, rowMap, colMap, NumEntriesPerRow),
151  myFirstRow_(0),
152  myNumRows_(0),
153  ignoreNonLocalEntries_(ignoreNonLocalEntries),
154 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
157 #endif
158 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
161 #endif
162  nonlocalCoefs_(),
163  workData_(128),
164  useNonlocalMatrix_ (false),
165  nonlocalMatrix_ (NULL),
166  sourceMap_(NULL),
167  colMap_(NULL),
168  exporter_(NULL),
169  tempMat_(NULL)
170 {
171  myFirstRow_ = rowMap.MinMyGID64();
172  myNumRows_ = rowMap.NumMyElements();
173 }
174 
175 //----------------------------------------------------------------------------
177  const Epetra_CrsGraph& graph,
178  bool ignoreNonLocalEntries)
179  : Epetra_CrsMatrix(CV, graph),
180  myFirstRow_(0),
181  myNumRows_(0),
182  ignoreNonLocalEntries_(ignoreNonLocalEntries),
183 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
186 #endif
187 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
190 #endif
191  nonlocalCoefs_(),
192  workData_(128),
193  useNonlocalMatrix_ (false),
194  nonlocalMatrix_ (NULL),
195  sourceMap_(NULL),
196  colMap_(NULL),
197  exporter_(NULL),
198  tempMat_(NULL)
199 {
202 }
203 
204 //----------------------------------------------------------------------------
206  const Epetra_FECrsGraph& graph,
207  bool ignoreNonLocalEntries)
208  : Epetra_CrsMatrix(CV, graph),
209  myFirstRow_(0),
210  myNumRows_(0),
211  ignoreNonLocalEntries_(ignoreNonLocalEntries),
212 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
215 #endif
216 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
219 #endif
220  nonlocalCoefs_(),
221  workData_(128),
222  useNonlocalMatrix_ (graph.UseNonlocalGraph() && graph.nonlocalGraph_ != 0),
224  new Epetra_CrsMatrix(Copy,*graph.nonlocalGraph_) : NULL),
225  sourceMap_(NULL),
226  colMap_(NULL),
227  exporter_(NULL),
228  tempMat_(NULL)
229 {
232 }
233 
234 //----------------------------------------------------------------------------
236  : Epetra_CrsMatrix(src),
237  myFirstRow_(0),
238  myNumRows_(0),
239  ignoreNonLocalEntries_(false),
240 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
243 #endif
244 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
247 #endif
248  nonlocalCoefs_(),
249  workData_(128),
250  nonlocalMatrix_ (NULL),
251  sourceMap_(NULL),
252  colMap_(NULL),
253  exporter_(NULL),
254  tempMat_(NULL)
255 {
256  operator=(src);
257 }
258 
259 //----------------------------------------------------------------------------
261 {
262  if (this == &src) {
263  return( *this );
264  }
265 
266  DeleteMemory();
267 
269 
271 
272  myFirstRow_ = src.myFirstRow_;
273  myNumRows_ = src.myNumRows_;
275 
276 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
277  if (src.RowMap().GlobalIndicesInt() && nonlocalRows_int_.size() < 1) {
278  return( *this );
279  }
280 #endif
281 
282 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
283  if (src.RowMap().GlobalIndicesLongLong() && nonlocalRows_LL_.size() < 1) {
284  return( *this );
285  }
286 #endif
287 
288  if (useNonlocalMatrix_ && src.nonlocalMatrix_ != 0) {
290  return( *this );
291  }
292 
293 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
294  if (src.RowMap().GlobalIndicesInt()) {
297  }
298 #endif
299 
300 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
301  if (src.RowMap().GlobalIndicesLongLong()) {
304  }
305 #endif
306 
308 
309  return( *this );
310 }
311 
312 //----------------------------------------------------------------------------
314 {
315  DeleteMemory();
316 }
317 
318 void Epetra_FECrsMatrix::Print(std::ostream& os) const
319 {
321 
322  if (ignoreNonLocalEntries_ || RowMap().Comm().NumProc()==1) return;
323 
324  int MyPID = RowMap().Comm().MyPID();
325  int NumProc = RowMap().Comm().NumProc();
326 
327  if (useNonlocalMatrix_) {
328  if (MyPID==0) {
329  os << "[FECrsMatrix] Nonlocal matrix:\n";
330  }
331  nonlocalMatrix_->Print(os);
332  return;
333  }
334 
335  if (MyPID==0) {
336  os << "[FECrsMatrix] Nonlocal rows:\n";
337  os.width(8);
338  os << " Processor ";
339  os.width(10);
340  os << " Row Index ";
341  os.width(10);
342  os << " Col Index ";
343  os.width(20);
344  os << " Value ";
345  os << std::endl;
346  }
347  for (int iproc=0; iproc < NumProc; iproc++) {
348  if (MyPID==iproc) {
349  if(RowMap().GlobalIndicesInt()) {
350  const int nnr = nonlocalRows_int_.size();
351  for (int i=0; i<nnr; ++i) {
352  const int Row = nonlocalRows_int_[i];
353  const int ncols = nonlocalCols_int_[i].size();;
354  for (int j=0; j<ncols; ++j) {
355  os.width(8);
356  os << MyPID ; os << " ";
357  os.width(10);
358  os << Row ; os << " ";
359  os.width(10);
360  os << nonlocalCols_int_[i][j]; os << " ";
361  os.width(20);
362  os << nonlocalCoefs_[i][j]; os << " ";
363  os << std::endl;
364  }
365  }
366  } else {
367  const int nnr = nonlocalRows_LL_.size();
368  for (int i=0; i<nnr; ++i) {
369  const int Row = nonlocalRows_LL_[i];
370  const int ncols = nonlocalCols_LL_[i].size();;
371  for (int j=0; j<ncols; ++j) {
372  os.width(8);
373  os << MyPID ; os << " ";
374  os.width(10);
375  os << Row ; os << " ";
376  os.width(10);
377  os << nonlocalCols_LL_[i][j]; os << " ";
378  os.width(20);
379  os << nonlocalCoefs_[i][j]; os << " ";
380  os << std::endl;
381  }
382  }
383  }
384  }
385  // Do a few global ops to give I/O a chance to complete
386  RowMap().Comm().Barrier();
387  RowMap().Comm().Barrier();
388  RowMap().Comm().Barrier();
389  }
390 }
391 
392 //----------------------------------------------------------------------------
394 {
395 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
396  if (RowMap().GlobalIndicesInt()) {
397  nonlocalRows_int_.clear();
398  nonlocalCols_int_.clear();
399  }
400 #endif
401 
402 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
403  if (RowMap().GlobalIndicesLongLong()) {
404  nonlocalRows_LL_.clear();
405  nonlocalCols_LL_.clear();
406  }
407 #endif
408 
409  nonlocalCoefs_.clear();
410 
411  if (nonlocalMatrix_ != 0)
412  delete nonlocalMatrix_;
413 
414  if ( sourceMap_ )
415  delete sourceMap_;
416  if ( colMap_ )
417  delete colMap_;
418  if ( exporter_ )
419  delete exporter_;
420  if ( tempMat_ )
421  delete tempMat_;
422 
423 }
424 
425 //----------------------------------------------------------------------------
426 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
427 int Epetra_FECrsMatrix::SumIntoGlobalValues(int numIndices, const int* indices,
428  const double* const* values,
429  int format)
430 {
431  return(InputGlobalValues(numIndices, indices,
432  numIndices, indices,
433  values, format, SUMINTO));
434 }
435 #endif
436 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
437 int Epetra_FECrsMatrix::SumIntoGlobalValues(int numIndices, const long long* indices,
438  const double* const* values,
439  int format)
440 {
441  return(InputGlobalValues(numIndices, indices,
442  numIndices, indices,
443  values, format, SUMINTO));
444 }
445 #endif
446 //----------------------------------------------------------------------------
447 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
448 int Epetra_FECrsMatrix::SumIntoGlobalValues(int numRows, const int* rows,
449  int numCols, const int* cols,
450  const double* const* values,
451  int format)
452 {
453  return(InputGlobalValues(numRows, rows,
454  numCols, cols,
455  values, format, SUMINTO));
456 }
457 #endif
458 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
459 int Epetra_FECrsMatrix::SumIntoGlobalValues(int numRows, const long long* rows,
460  int numCols, const long long* cols,
461  const double* const* values,
462  int format)
463 {
464  return(InputGlobalValues(numRows, rows,
465  numCols, cols,
466  values, format, SUMINTO));
467 }
468 #endif
469 //----------------------------------------------------------------------------
470 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
471 int Epetra_FECrsMatrix::SumIntoGlobalValues(int numIndices, const int* indices,
472  const double* values,
473  int format)
474 {
475  return(InputGlobalValues(numIndices, indices,
476  numIndices, indices,
477  values, format, SUMINTO));
478 }
479 #endif
480 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
481 int Epetra_FECrsMatrix::SumIntoGlobalValues(int numIndices, const long long* indices,
482  const double* values,
483  int format)
484 {
485  return(InputGlobalValues(numIndices, indices,
486  numIndices, indices,
487  values, format, SUMINTO));
488 }
489 #endif
490 //----------------------------------------------------------------------------
491 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
492 int Epetra_FECrsMatrix::SumIntoGlobalValues(int numRows, const int* rows,
493  int numCols, const int* cols,
494  const double* values,
495  int format)
496 {
497  return(InputGlobalValues(numRows, rows,
498  numCols, cols,
499  values, format, SUMINTO));
500 }
501 
502 #endif
503 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
504 int Epetra_FECrsMatrix::SumIntoGlobalValues(int numRows, const long long* rows,
505  int numCols, const long long* cols,
506  const double* values,
507  int format)
508 {
509  return(InputGlobalValues(numRows, rows,
510  numCols, cols,
511  values, format, SUMINTO));
512 }
513 #endif
514 //----------------------------------------------------------------------------
515 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
517  const Epetra_SerialDenseMatrix& values,
518  int format)
519 {
520  if (indices.Length() != values.M() || indices.Length() != values.N()) {
521  return(-1);
522  }
523 
524  return( SumIntoGlobalValues(indices.Length(), indices.Values(),
525  values.A(), format) );
526 }
527 
528 #endif
529 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
531  const Epetra_SerialDenseMatrix& values,
532  int format)
533 {
534  if (indices.Length() != values.M() || indices.Length() != values.N()) {
535  return(-1);
536  }
537 
538  return( SumIntoGlobalValues(indices.Length(), indices.Values(),
539  values.A(), format) );
540 }
541 #endif
542 //----------------------------------------------------------------------------
543 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
545  const Epetra_SerialDenseMatrix& values,
546  int format)
547 {
548  if (indices.Length() != values.M() || indices.Length() != values.N()) {
549  return(-1);
550  }
551 
552  return( InsertGlobalValues(indices.Length(), indices.Values(),
553  values.A(), format) );
554 }
555 #endif
556 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
558  const Epetra_SerialDenseMatrix& values,
559  int format)
560 {
561  if (indices.Length() != values.M() || indices.Length() != values.N()) {
562  return(-1);
563  }
564 
565  return( InsertGlobalValues(indices.Length(), indices.Values(),
566  values.A(), format) );
567 }
568 #endif
569 //----------------------------------------------------------------------------
570 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
572  const Epetra_SerialDenseMatrix& values,
573  int format)
574 {
575  if (indices.Length() != values.M() || indices.Length() != values.N()) {
576  return(-1);
577  }
578 
579  return( ReplaceGlobalValues(indices.Length(), indices.Values(),
580  values.A(), format) );
581 }
582 #endif
583 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
585  const Epetra_SerialDenseMatrix& values,
586  int format)
587 {
588  if (indices.Length() != values.M() || indices.Length() != values.N()) {
589  return(-1);
590  }
591 
592  return( ReplaceGlobalValues(indices.Length(), indices.Values(),
593  values.A(), format) );
594 }
595 #endif
596 //----------------------------------------------------------------------------
597 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
599  const Epetra_IntSerialDenseVector& cols,
600  const Epetra_SerialDenseMatrix& values,
601  int format)
602 {
603  if (rows.Length() != values.M() || cols.Length() != values.N()) {
604  return(-1);
605  }
606 
607  return( SumIntoGlobalValues(rows.Length(), rows.Values(),
608  cols.Length(), cols.Values(),
609  values.A(), format) );
610 }
611 #endif
612 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
615  const Epetra_SerialDenseMatrix& values,
616  int format)
617 {
618  if (rows.Length() != values.M() || cols.Length() != values.N()) {
619  return(-1);
620  }
621 
622  return( SumIntoGlobalValues(rows.Length(), rows.Values(),
623  cols.Length(), cols.Values(),
624  values.A(), format) );
625 }
626 #endif
627 //----------------------------------------------------------------------------
628 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
630  const Epetra_IntSerialDenseVector& cols,
631  const Epetra_SerialDenseMatrix& values,
632  int format)
633 {
634  if (rows.Length() != values.M() || cols.Length() != values.N()) {
635  return(-1);
636  }
637 
638  return( InsertGlobalValues(rows.Length(), rows.Values(),
639  cols.Length(), cols.Values(),
640  values.A(), format) );
641 }
642 #endif
643 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
646  const Epetra_SerialDenseMatrix& values,
647  int format)
648 {
649  if (rows.Length() != values.M() || cols.Length() != values.N()) {
650  return(-1);
651  }
652 
653  return( InsertGlobalValues(rows.Length(), rows.Values(),
654  cols.Length(), cols.Values(),
655  values.A(), format) );
656 }
657 #endif
658 //----------------------------------------------------------------------------
659 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
661  const Epetra_IntSerialDenseVector& cols,
662  const Epetra_SerialDenseMatrix& values,
663  int format)
664 {
665  if (rows.Length() != values.M() || cols.Length() != values.N()) {
666  return(-1);
667  }
668 
669  return( ReplaceGlobalValues(rows.Length(), rows.Values(),
670  cols.Length(), cols.Values(),
671  values.A(), format) );
672 }
673 #endif
674 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
675 
678  const Epetra_SerialDenseMatrix& values,
679  int format)
680 {
681  if (rows.Length() != values.M() || cols.Length() != values.N()) {
682  return(-1);
683  }
684 
685  return( ReplaceGlobalValues(rows.Length(), rows.Values(),
686  cols.Length(), cols.Values(),
687  values.A(), format) );
688 }
689 #endif
690 //----------------------------------------------------------------------------
691 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
692 int Epetra_FECrsMatrix::InsertGlobalValues(int numIndices, const int* indices,
693  const double* const* values,
694  int format)
695 {
696  return(InputGlobalValues(numIndices, indices,
697  numIndices, indices,
698  values, format, INSERT));
699 }
700 #endif
701 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
702 
703 int Epetra_FECrsMatrix::InsertGlobalValues(int numIndices, const long long* indices,
704  const double* const* values,
705  int format)
706 {
707  return(InputGlobalValues(numIndices, indices,
708  numIndices, indices,
709  values, format, INSERT));
710 }
711 #endif
712 //----------------------------------------------------------------------------
713 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
714 int Epetra_FECrsMatrix::InsertGlobalValues(int numRows, const int* rows,
715  int numCols, const int* cols,
716  const double* const* values,
717  int format)
718 {
719  return(InputGlobalValues(numRows, rows,
720  numCols, cols,
721  values, format, INSERT));
722 }
723 #endif
724 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
725 
726 int Epetra_FECrsMatrix::InsertGlobalValues(int numRows, const long long* rows,
727  int numCols, const long long* cols,
728  const double* const* values,
729  int format)
730 {
731  return(InputGlobalValues(numRows, rows,
732  numCols, cols,
733  values, format, INSERT));
734 }
735 #endif
736 //----------------------------------------------------------------------------
737 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
738 int Epetra_FECrsMatrix::InsertGlobalValues(int numIndices, const int* indices,
739  const double* values,
740  int format)
741 {
742  return(InputGlobalValues(numIndices, indices,
743  numIndices, indices,
744  values, format, INSERT));
745 }
746 #endif
747 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
748 
749 int Epetra_FECrsMatrix::InsertGlobalValues(int numIndices, const long long* indices,
750  const double* values,
751  int format)
752 {
753  return(InputGlobalValues(numIndices, indices,
754  numIndices, indices,
755  values, format, INSERT));
756 }
757 #endif
758 //----------------------------------------------------------------------------
759 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
760 int Epetra_FECrsMatrix::InsertGlobalValues(int numRows, const int* rows,
761  int numCols, const int* cols,
762  const double* values,
763  int format)
764 {
765  return(InputGlobalValues(numRows, rows,
766  numCols, cols,
767  values, format, INSERT));
768 }
769 #endif
770 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
771 
772 int Epetra_FECrsMatrix::InsertGlobalValues(int numRows, const long long* rows,
773  int numCols, const long long* cols,
774  const double* values,
775  int format)
776 {
777  return(InputGlobalValues(numRows, rows,
778  numCols, cols,
779  values, format, INSERT));
780 }
781 #endif
782 //----------------------------------------------------------------------------
783 template<typename int_type>
784 int Epetra_FECrsMatrix::SumIntoGlobalValues(int_type GlobalRow, int NumEntries,
785  const double* values, const int_type* Indices)
786 {
787  if (Map().MyGID(GlobalRow))
788  return Epetra_CrsMatrix::SumIntoGlobalValues(GlobalRow, NumEntries,
789  values, Indices);
790  else if (useNonlocalMatrix_)
791  return nonlocalMatrix_->SumIntoGlobalValues(GlobalRow,
792  NumEntries, values, Indices);
793  else
794  return InputNonlocalGlobalValues(GlobalRow, NumEntries, Indices, values, SUMINTO);
795 }
796 
797 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
798 int Epetra_FECrsMatrix::SumIntoGlobalValues(int GlobalRow, int NumEntries,
799  const double* values, const int* Indices)
800 {
801  if(RowMap().GlobalIndicesInt())
802  return SumIntoGlobalValues<int>(GlobalRow, NumEntries, values, Indices);
803  else
804  throw ReportError("Epetra_FECrsMatrix::SumIntoGlobalValues int version called for a matrix that is not int.", -1);
805 }
806 #endif
807 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
808 
809 
810 int Epetra_FECrsMatrix::SumIntoGlobalValues(long long GlobalRow, int NumEntries,
811  const double* values, const long long* Indices)
812 {
813  if(RowMap().GlobalIndicesLongLong())
814  return SumIntoGlobalValues<long long>(GlobalRow, NumEntries, values, Indices);
815  else
816  throw ReportError("Epetra_FECrsMatrix::SumIntoGlobalValues long long version called for a matrix that is not long long.", -1);
817 }
818 #endif
819 //----------------------------------------------------------------------------
820 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
821 int Epetra_FECrsMatrix::InsertGlobalValues(int GlobalRow, int NumEntries,
822  const double* values, const int* Indices)
823 {
824  return(InputGlobalValues(1, &GlobalRow,
825  NumEntries, Indices, values,
827 }
828 #endif
829 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
830 
831 int Epetra_FECrsMatrix::InsertGlobalValues(long long GlobalRow, int NumEntries,
832  const double* values, const long long* Indices)
833 {
834  return(InputGlobalValues(1, &GlobalRow,
835  NumEntries, Indices, values,
837 }
838 #endif
839 //----------------------------------------------------------------------------
840 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
841 int Epetra_FECrsMatrix::InsertGlobalValues(int GlobalRow, int NumEntries,
842  double* values, int* Indices)
843 {
844  return(InputGlobalValues(1, &GlobalRow,
845  NumEntries, Indices, values,
847 }
848 #endif
849 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
850 int Epetra_FECrsMatrix::InsertGlobalValues(long long GlobalRow, int NumEntries,
851  double* values, long long* Indices)
852 {
853  return(InputGlobalValues(1, &GlobalRow,
854  NumEntries, Indices, values,
856 }
857 #endif
858 //----------------------------------------------------------------------------
859 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
860 int Epetra_FECrsMatrix::ReplaceGlobalValues(int GlobalRow, int NumEntries,
861  const double* values, const int* Indices)
862 {
863  return(InputGlobalValues(1, &GlobalRow,
864  NumEntries, Indices, values,
866 }
867 #endif
868 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
869 int Epetra_FECrsMatrix::ReplaceGlobalValues(long long GlobalRow, int NumEntries,
870  const double* values, const long long* Indices)
871 {
872  return(InputGlobalValues(1, &GlobalRow,
873  NumEntries, Indices, values,
875 }
876 #endif
877 //----------------------------------------------------------------------------
878 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
879 int Epetra_FECrsMatrix::ReplaceGlobalValues(int numIndices, const int* indices,
880  const double* const* values,
881  int format)
882 {
883  return(InputGlobalValues(numIndices, indices,
884  numIndices, indices,
885  values, format, REPLACE));
886 }
887 #endif
888 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
889 int Epetra_FECrsMatrix::ReplaceGlobalValues(int numIndices, const long long* indices,
890  const double* const* values,
891  int format)
892 {
893  return(InputGlobalValues(numIndices, indices,
894  numIndices, indices,
895  values, format, REPLACE));
896 }
897 #endif
898 //----------------------------------------------------------------------------
899 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
900 int Epetra_FECrsMatrix::ReplaceGlobalValues(int numRows, const int* rows,
901  int numCols, const int* cols,
902  const double* const* values,
903  int format)
904 {
905  return(InputGlobalValues(numRows, rows,
906  numCols, cols,
907  values, format, REPLACE));
908 }
909 #endif
910 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
911 int Epetra_FECrsMatrix::ReplaceGlobalValues(int numRows, const long long* rows,
912  int numCols, const long long* cols,
913  const double* const* values,
914  int format)
915 {
916  return(InputGlobalValues(numRows, rows,
917  numCols, cols,
918  values, format, REPLACE));
919 }
920 #endif
921 //----------------------------------------------------------------------------
922 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
923 int Epetra_FECrsMatrix::ReplaceGlobalValues(int numIndices, const int* indices,
924  const double* values,
925  int format)
926 {
927  return(InputGlobalValues(numIndices, indices,
928  numIndices, indices,
929  values, format, REPLACE));
930 }
931 #endif
932 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
933 int Epetra_FECrsMatrix::ReplaceGlobalValues(int numIndices, const long long* indices,
934  const double* values,
935  int format)
936 {
937  return(InputGlobalValues(numIndices, indices,
938  numIndices, indices,
939  values, format, REPLACE));
940 }
941 #endif
942 //----------------------------------------------------------------------------
943 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
944 int Epetra_FECrsMatrix::ReplaceGlobalValues(int numRows, const int* rows,
945  int numCols, const int* cols,
946  const double* values,
947  int format)
948 {
949  return(InputGlobalValues(numRows, rows,
950  numCols, cols,
951  values, format, REPLACE));
952 }
953 #endif
954 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
955 int Epetra_FECrsMatrix::ReplaceGlobalValues(int numRows, const long long* rows,
956  int numCols, const long long* cols,
957  const double* values,
958  int format)
959 {
960  return(InputGlobalValues(numRows, rows,
961  numCols, cols,
962  values, format, REPLACE));
963 }
964 #endif
965 //----------------------------------------------------------------------------
966 int Epetra_FECrsMatrix::GlobalAssemble(bool callFillComplete, Epetra_CombineMode combineMode,
967  bool save_off_and_reuse_map_exporter)
968 {
969  return( GlobalAssemble(DomainMap(), RangeMap(), callFillComplete, combineMode, save_off_and_reuse_map_exporter));
970 }
971 
972 //----------------------------------------------------------------------------
973 template<typename int_type>
975  const Epetra_Map& range_map,
976  bool callFillComplete,
977  Epetra_CombineMode combineMode,
978  bool save_off_and_reuse_map_exporter)
979 {
980  if (Map().Comm().NumProc() < 2 || ignoreNonLocalEntries_) {
981  if (callFillComplete) {
982  EPETRA_CHK_ERR( FillComplete(domain_map, range_map) );
983  }
984  return(0);
985  }
986 
987  std::vector<int_type>& nonlocalRows_var = nonlocalRows<int_type>();
988  std::vector<std::vector<int_type> >& nonlocalCols_var = nonlocalCols<int_type>();
989 
990  if (useNonlocalMatrix_) {
992  }
993  else {
994  //In this method we need to gather all the non-local (overlapping) data
995  //that's been input on each processor, into the
996  //non-overlapping distribution defined by the map that 'this' matrix was
997  //constructed with.
998 
999  //First build a map that describes our nonlocal data.
1000  //We'll use the arbitrary distribution constructor of Map.
1001 
1002  int_type* nlr_ptr = Epetra_Util_data_ptr(nonlocalRows_var);
1003  if (sourceMap_ == NULL)
1004  sourceMap_ = new Epetra_Map((int_type) -1, (int) nonlocalRows_var.size(), nlr_ptr,
1005  (int_type) Map().IndexBase64(), Map().Comm());
1006 
1007  //If sourceMap has global size 0, then no nonlocal data exists and we can
1008  //skip most of this function.
1009  if (sourceMap_->NumGlobalElements64() < 1) {
1010  if (callFillComplete) {
1011  EPETRA_CHK_ERR( FillComplete(domain_map, range_map) );
1012  }
1013  if (!save_off_and_reuse_map_exporter) {
1014  delete sourceMap_;
1015  sourceMap_ = NULL;
1016  }
1017  return(0);
1018  }
1019 
1020  //We also need to build a column-map, containing the columns in our
1021  //nonlocal data. To do that, create a list of all column-indices that
1022  //occur in our nonlocal rows.
1023  bool first_time=!save_off_and_reuse_map_exporter;
1024  if ( colMap_ == NULL ) {
1025  first_time = true;
1026  std::vector<int_type> cols;
1027 
1028  for(size_t i=0; i<nonlocalRows_var.size(); ++i) {
1029  for(size_t j=0; j<nonlocalCols_var[i].size(); ++j) {
1030  int_type col = nonlocalCols_var[i][j];
1031  typename std::vector<int_type>::iterator it =
1032  std::lower_bound(cols.begin(), cols.end(), col);
1033  if (it == cols.end() || *it != col) {
1034  cols.insert(it, col);
1035  }
1036  }
1037  }
1038 
1039  int_type* cols_ptr = Epetra_Util_data_ptr(cols);
1040 
1041  colMap_ = new Epetra_Map((int_type) -1, (int) cols.size(), cols_ptr,
1042  (int_type) Map().IndexBase64(), Map().Comm());
1043  }
1044  //now we need to create a matrix with sourceMap and colMap, and fill it with
1045  //our nonlocal data so we can then export it to the correct owning processors.
1046 
1047  std::vector<int> nonlocalRowLengths(nonlocalRows_var.size());
1048  for(size_t i=0; i<nonlocalRows_var.size(); ++i) {
1049  nonlocalRowLengths[i] = (int) nonlocalCols_var[i].size();
1050  }
1051 
1052  int* nlRLptr = Epetra_Util_data_ptr(nonlocalRowLengths);
1053  if ( first_time && tempMat_ == NULL )
1054  tempMat_ = new Epetra_CrsMatrix(Copy, *sourceMap_, *colMap_, nlRLptr);
1055  else
1056  tempMat_->PutScalar(0.);
1057 
1058  for(size_t i=0; i<nonlocalRows_var.size(); ++i) {
1059  if ( first_time ) {
1060  EPETRA_CHK_ERR( tempMat_->InsertGlobalValues(nonlocalRows_var[i],
1061  (int) nonlocalCols_var[i].size(),
1063  Epetra_Util_data_ptr(nonlocalCols_var[i])) );
1064  } else {
1065  EPETRA_CHK_ERR( tempMat_->SumIntoGlobalValues(nonlocalRows_var[i],
1066  (int) nonlocalCols_var[i].size(),
1068  Epetra_Util_data_ptr(nonlocalCols_var[i])) );
1069  }
1070  }
1071 
1072  if (!save_off_and_reuse_map_exporter) {
1073  delete sourceMap_;
1074  delete colMap_;
1075  sourceMap_ = colMap_ = NULL;
1076  }
1077 
1078  //Next we need to make sure the 'indices-are-global' attribute of tempMat's
1079  //graph is set to true, in case this processor doesn't end up calling the
1080  //InsertGlobalValues method...
1081 
1082  if (first_time) {
1083  const Epetra_CrsGraph& graph = tempMat_->Graph();
1084  Epetra_CrsGraph& nonconst_graph = const_cast<Epetra_CrsGraph&>(graph);
1085  nonconst_graph.SetIndicesAreGlobal(true);
1086  }
1087  }
1088 
1089  //Now we need to call FillComplete on our temp matrix. We need to
1090  //pass a DomainMap and RangeMap, which are not the same as the RowMap
1091  //and ColMap that we constructed the matrix with.
1092  EPETRA_CHK_ERR(tempMat_->FillComplete(domain_map, range_map));
1093 
1094  if (exporter_ == NULL)
1096 
1097  EPETRA_CHK_ERR(Export(*tempMat_, *exporter_, combineMode));
1098 
1099  if(callFillComplete) {
1100  EPETRA_CHK_ERR(FillComplete(domain_map, range_map));
1101  }
1102 
1103  //now reset the values in our nonlocal data
1104  if (useNonlocalMatrix_) {
1105  nonlocalMatrix_->PutScalar(0.0);
1106  } else {
1107  for(size_t i=0; i<nonlocalRows_var.size(); ++i) {
1108  nonlocalCols_var[i].resize(0);
1109  nonlocalCoefs_[i].resize(0);
1110  }
1111  }
1112 
1113  if (!save_off_and_reuse_map_exporter) {
1114  delete exporter_;
1115  exporter_ = NULL;
1116  if (!useNonlocalMatrix_)
1117  delete tempMat_;
1118  tempMat_ = NULL;
1119  }
1120  return(0);
1121 }
1122 
1123 int Epetra_FECrsMatrix::GlobalAssemble(const Epetra_Map& domain_map,
1124  const Epetra_Map& range_map,
1125  bool callFillComplete,
1126  Epetra_CombineMode combineMode,
1127  bool save_off_and_reuse_map_exporter)
1128 {
1129  if(!domain_map.GlobalIndicesTypeMatch(range_map))
1130  throw ReportError("Epetra_FECrsMatrix::GlobalAssemble: cannot be called with different indices types for domainMap and rangeMap", -1);
1131 
1132  if(!RowMap().GlobalIndicesTypeMatch(domain_map))
1133  throw ReportError("Epetra_FECrsMatrix::GlobalAssemble: cannot be called with different indices types for row map and incoming rangeMap", -1);
1134 
1135  if(RowMap().GlobalIndicesInt())
1136 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
1137  return GlobalAssemble<int>(domain_map, range_map, callFillComplete, combineMode, save_off_and_reuse_map_exporter);
1138 #else
1139  throw ReportError("Epetra_FECrsMatrix::GlobalAssemble: ERROR, GlobalIndicesInt but no API for it.",-1);
1140 #endif
1141 
1142  if(RowMap().GlobalIndicesLongLong())
1143 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
1144  return GlobalAssemble<long long>(domain_map, range_map, callFillComplete, combineMode, save_off_and_reuse_map_exporter);
1145 #else
1146  throw ReportError("Epetra_FECrsMatrix::GlobalAssemble: ERROR, GlobalIndicesLongLong but no API for it.",-1);
1147 #endif
1148 
1149  throw ReportError("Epetra_FECrsMatrix::GlobalAssemble: Internal error, unable to determine global index type of maps", -1);
1150 }
1151 
1152 //----------------------------------------------------------------------------
1153 template<typename int_type>
1155  int numRows, const int_type* rows,
1156  int numCols, const int_type* cols,
1157  const double* values,
1158  int mode)
1159 {
1160  if(!RowMap().template GlobalIndicesIsType<int_type>())
1161  throw ReportError("Epetra_FECrsMatrix::InputGlobalValues_RowMajor mismatch between argument types (int/long long) and map type.", -1);
1162 
1163  int returncode = 0;
1164  int err = 0;
1165 
1166  for(int i=0; i<numRows; ++i) {
1167  double* valuesptr = (double*)values + i*numCols;
1168 
1169  int local_row_id = Map().LID(rows[i]);
1170  if (local_row_id >= 0) {
1171  switch(mode) {
1173  err = this->Epetra_CrsMatrix::SumIntoGlobalValues(rows[i], numCols,
1174  valuesptr, (int_type*)cols);
1175  if (err<0) return(err);
1176  if (err>0) returncode = err;
1177  break;
1179  err = this->Epetra_CrsMatrix::ReplaceGlobalValues(rows[i], numCols,
1180  valuesptr, (int_type*)cols);
1181  if (err<0) return(err);
1182  if (err>0) returncode = err;
1183  break;
1185  err = this->Epetra_CrsMatrix::InsertGlobalValues(rows[i], numCols,
1186  valuesptr, (int_type*)cols);
1187  if (err<0) return(err);
1188  if (err>0) returncode = err;
1189  break;
1190  default:
1191  std::cerr << "Epetra_FECrsMatrix: internal error, bad input mode."<< std::endl;
1192  return(-1);
1193  }
1194  }
1195  else {
1196 #ifdef EPETRA_HAVE_OMP
1197 #ifdef EPETRA_HAVE_OMP_NONASSOCIATIVE
1198  if (! ignoreNonLocalEntries_) {
1199 #endif
1200 #endif
1201  err = InputNonlocalGlobalValues(rows[i], numCols, cols,
1202  valuesptr, mode);
1203  if (err<0) return(err);
1204  if (err>0) returncode = err;
1205 #ifdef EPETRA_HAVE_OMP
1206 #ifdef EPETRA_HAVE_OMP_NONASSOCIATIVE
1207  }
1208 #endif
1209 #endif
1210  }
1211  }
1212 
1213  return(returncode);
1214 }
1215 
1216 //----------------------------------------------------------------------------
1217 template<typename int_type>
1218 int Epetra_FECrsMatrix::InputGlobalValues(int numRows, const int_type* rows,
1219  int numCols, const int_type* cols,
1220  const double*const* values,
1221  int format, int mode)
1222 {
1223  if(!RowMap().template GlobalIndicesIsType<int_type>())
1224  throw ReportError("Epetra_FECrsMatrix::InputGlobalValues mismatch between argument types (int/long long) and map type.", -1);
1225 
1226  if (format != Epetra_FECrsMatrix::ROW_MAJOR &&
1228  std::cerr << "Epetra_FECrsMatrix: unrecognized format specifier."<< std::endl;
1229  return(-1);
1230  }
1231 
1232  if (format == Epetra_FECrsMatrix::COLUMN_MAJOR) {
1233  workData_.resize(numCols);
1234  }
1235 
1236  int returncode = 0;
1237 
1238  for(int i=0; i<numRows; ++i) {
1239  if (format == Epetra_FECrsMatrix::ROW_MAJOR) {
1240  returncode += InputGlobalValues_RowMajor(1, &rows[i], numCols, cols,
1241  values[i], mode);
1242  if (returncode < 0) return returncode;
1243  continue;
1244  }
1245 
1246  //If we get to here, the data is in column-major order.
1247 
1248  double* valuesptr = Epetra_Util_data_ptr(workData_);
1249 
1250  //Since the data is in column-major order, then we copy the i-th row
1251  //of the values table into workData_, in order to have the row in
1252  //contiguous memory.
1253  //This is slow and not thread-safe.
1254 
1255  for(int j=0; j<numCols; ++j) {
1256  valuesptr[j] = values[j][i];
1257  }
1258 
1259  returncode += InputGlobalValues_RowMajor(1, &rows[i], numCols, cols, valuesptr, mode);
1260  if (returncode < 0) return returncode;
1261  }
1262 
1263  return(returncode);
1264 }
1265 
1266 //----------------------------------------------------------------------------
1267 template<typename int_type>
1268 int Epetra_FECrsMatrix::InputGlobalValues(int numRows, const int_type* rows,
1269  int numCols, const int_type* cols,
1270  const double* values,
1271  int format, int mode)
1272 {
1273  if(!RowMap().template GlobalIndicesIsType<int_type>())
1274  throw ReportError("Epetra_FECrsMatrix::InputGlobalValues mismatch between argument types (int/long long) and map type.", -1);
1275 
1276  if (format == Epetra_FECrsMatrix::ROW_MAJOR) {
1277  return InputGlobalValues_RowMajor(numRows, rows, numCols, cols, values, mode);
1278  }
1279 
1280  workData_.resize(numCols);
1281 
1282  int returncode = 0;
1283  for(int i=0; i<numRows; ++i) {
1284  //copy each row out of the column-major values array, so we can pass it
1285  //to a row-major input function.
1286  for(int j=0; j<numCols; ++j) {
1287  workData_[j] = values[i+j*numRows];
1288  }
1289  int err = InputGlobalValues_RowMajor(1, &rows[i], numCols, cols, Epetra_Util_data_ptr(workData_), mode);
1290  if (err < 0) return err;
1291  returncode += err;
1292  }
1293 
1294  return(returncode);
1295 }
1296 
1297 //----------------------------------------------------------------------------
1298 template<typename int_type>
1300  int numCols, const int_type* cols,
1301  const double* values,
1302  int mode)
1303 {
1304  if(!RowMap().template GlobalIndicesIsType<int_type>())
1305  throw ReportError("Epetra_FECrsMatrix::InputNonlocalGlobalValues mismatch between argument types (int/long long) and map type.", -1);
1306 
1307  // if we already have a nonlocal matrix object, this is easier...
1308  if (useNonlocalMatrix_) {
1309  int err, returncode = 0;
1310  double* valuesptr = (double*)values;
1311  switch(mode) {
1313  err = nonlocalMatrix_->SumIntoGlobalValues(row, numCols,
1314  valuesptr, (int_type*)cols);
1315  if (err<0) return(err);
1316  if (err>0) returncode = err;
1317  break;
1319  err = nonlocalMatrix_->ReplaceGlobalValues(row, numCols,
1320  valuesptr, (int_type*)cols);
1321  if (err<0) return(err);
1322  if (err>0) returncode = err;
1323  break;
1325  err = nonlocalMatrix_->InsertGlobalValues(row, numCols,
1326  valuesptr, (int_type*)cols);
1327  if (err<0) return(err);
1328  if (err>0) returncode = err;
1329  break;
1330  default:
1331  std::cerr << "Epetra_FECrsMatrix: internal error, bad input mode."<< std::endl;
1332  return(-1);
1333  }
1334  return (returncode);
1335  }
1336  int ierr1 = 0, ierr2 = 0;
1337 #ifdef EPETRA_HAVE_OMP
1338 #ifdef EPETRA_HAVE_OMP_NONASSOCIATIVE
1339 #pragma omp critical
1340 #endif
1341 #endif
1342  {
1343  std::vector<int_type>& nonlocalRows_var = nonlocalRows<int_type>();
1344 
1345  //find offset of this row in our list of nonlocal rows.
1346  typename std::vector<int_type>::iterator it =
1347  std::lower_bound(nonlocalRows_var.begin(), nonlocalRows_var.end(), row);
1348 
1349  int rowoffset = (int) (it - nonlocalRows_var.begin());
1350  if (it == nonlocalRows_var.end() || *it != row) {
1351  ierr1 = InsertNonlocalRow(row, it);
1352  }
1353 
1354  for(int i=0; i<numCols; ++i) {
1355  ierr2 = InputNonlocalValue(rowoffset, cols[i], values[i], mode);
1356  }
1357  }
1358  EPETRA_CHK_ERR(ierr1);
1359  EPETRA_CHK_ERR(ierr2);
1360 
1361  return(0);
1362 }
1363 
1364 //----------------------------------------------------------------------------
1365 template<typename int_type>
1366 int Epetra_FECrsMatrix::InsertNonlocalRow(int_type row, typename std::vector<int_type>::iterator iter)
1367 {
1368  if(!RowMap().template GlobalIndicesIsType<int_type>())
1369  throw ReportError("Epetra_FECrsMatrix::InsertNonlocalRow mismatch between argument types (int/long long) and map type.", -1);
1370 
1371  std::vector<int_type>& nonlocalRows_var = nonlocalRows<int_type>();
1372  std::vector<std::vector<int_type> >& nonlocalCols_var = nonlocalCols<int_type>();
1373 
1374  int offset = (int) (iter - nonlocalRows_var.begin());
1375  nonlocalRows_var.insert(iter, row);
1376  typename std::vector<std::vector<int_type> >::iterator cols_iter = nonlocalCols_var.begin() + offset;
1377  nonlocalCols_var.insert(cols_iter, std::vector<int_type>());
1378  std::vector<std::vector<double> >::iterator coefs_iter = nonlocalCoefs_.begin() + offset;
1379  nonlocalCoefs_.insert(coefs_iter, std::vector<double>());
1380 
1381  return(0);
1382 }
1383 
1384 //----------------------------------------------------------------------------
1385 template<typename int_type>
1387  int_type col, double value,
1388  int mode)
1389 {
1390  if(!RowMap().template GlobalIndicesIsType<int_type>())
1391  throw ReportError("Epetra_FECrsMatrix::InputNonlocalValue mismatch between argument types (int/long long) and map type.", -1);
1392 
1393  std::vector<int_type>& colIndices = nonlocalCols<int_type>()[rowoffset];
1394  std::vector<double>& coefs = nonlocalCoefs_[rowoffset];
1395 
1396  typename std::vector<int_type>::iterator it =
1397  std::lower_bound(colIndices.begin(), colIndices.end(), col);
1398 
1399  if (it == colIndices.end() || *it != col) {
1400  int offset = (int) (it - colIndices.begin());
1401  colIndices.insert(it, col);
1402  std::vector<double>::iterator dit = coefs.begin()+offset;
1403  coefs.insert(dit, value);
1404  return 0;
1405  }
1406 
1407  int coloffset = (int) (it - colIndices.begin());
1408  if (mode == SUMINTO || mode == INSERT) {
1409  coefs[coloffset] += value;
1410  }
1411  else {
1412  coefs[coloffset] = value;
1413  }
1414 
1415  return(0);
1416 }
1417 
long long MinMyGID64() const
Epetra_Map: A class for partitioning vectors and matrices.
Definition: Epetra_Map.h:119
const Epetra_Map & RangeMap() const
Returns the Epetra_Map object associated with the range of this matrix operator.
Epetra_IntSerialDenseVector: A class for constructing and using dense vectors.
int Length() const
Returns length of vector.
int InputGlobalValues_RowMajor(int numRows, const int_type *rows, int numCols, const int_type *cols, const double *values, int mode)
bool GlobalIndicesLongLong() const
Returns true if map create with long long NumGlobalElements.
long long * Values()
Returns pointer to the values in vector.
long long IndexBase64() const
Epetra_CrsMatrix * tempMat_
long long NumGlobalElements64() const
std::vector< std::vector< int > > nonlocalCols_int_
int InsertNonlocalRow(int_type row, typename std::vector< int_type >::iterator offset)
virtual int InsertGlobalValues(int GlobalRow, int NumEntries, const double *Values, const int *Indices)
Insert a list of elements in a given global row of the matrix.
Epetra Finite-Element CrsGraph.
virtual int SumIntoGlobalValues(int GlobalRow, int NumEntries, const double *Values, const int *Indices)
Add this list of entries to existing values for a given global row of the matrix. ...
T * Epetra_Util_data_ptr(std::vector< T > &vec)
Function that returns either a pointer to the first entry in the vector or, if the vector is empty...
Definition: Epetra_Util.h:422
#define EPETRA_CHK_ERR(a)
Epetra_Export: This class builds an export object for efficient exporting of off-processor elements...
Definition: Epetra_Export.h:62
int InputNonlocalGlobalValues(int_type row, int numCols, const int_type *cols, const double *values, int mode)
bool GlobalIndicesInt() const
Returns true if map create with int NumGlobalElements.
Epetra_SerialDenseMatrix: A class for constructing and using real double precision general dense matr...
virtual void Barrier() const =0
Epetra_Comm Barrier function.
std::vector< int > nonlocalRows_int_
virtual int MyPID() const =0
Return my process ID.
int FillComplete(bool OptimizeDataStorage=true)
Signal that data entry is complete. Perform transformations to local index space. ...
int Length() const
Returns length of vector.
int InsertGlobalValues(int GlobalRow, int NumEntries, const double *Values, const int *Indices)
override base-class Epetra_CrsMatrix::InsertGlobalValues method
std::vector< long long > nonlocalRows_LL_
int PutScalar(double ScalarConstant)
Initialize all values in the matrix with constant value.
double * A() const
Returns pointer to the this matrix.
const Epetra_Map & RowMap() const
Returns the Epetra_Map object associated with the rows of this matrix.
int NumMyElements() const
Number of elements on the calling processor.
void Print(std::ostream &os) const
Print method.
Epetra Finite-Element CrsMatrix.
Epetra_FECrsMatrix & operator=(const Epetra_FECrsMatrix &src)
Assignment operator.
int Export(const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
Exports an Epetra_DistObject using the Epetra_Import object.
void SetIndicesAreGlobal(bool Flag)
Epetra_Export * exporter_
int GlobalAssemble(bool callFillComplete=true, Epetra_CombineMode combineMode=Add, bool save_off_and_reuse_map_exporter=false)
Gather any overlapping/shared data into the non-overlapping partitioning defined by the Map that was ...
Epetra_LongLongSerialDenseVector: A class for constructing and using dense vectors.
int InputGlobalValues(int numRows, const int_type *rows, int numCols, const int_type *cols, const double *const *values, int format, int mode)
virtual void Print(std::ostream &os) const
Print method.
int LID(int GID) const
Returns local ID of global ID, return -1 if not found on this processor.
std::vector< std::vector< long long > > nonlocalCols_LL_
const Epetra_Comm & Comm() const
Access function for Epetra_Comm communicator.
std::vector< std::vector< double > > nonlocalCoefs_
int SumIntoGlobalValues(int GlobalRow, int NumEntries, const double *Values, const int *Indices)
override base-class Epetra_CrsMatrix::SumIntoGlobalValues method
virtual int NumProc() const =0
Returns total number of processes.
Epetra_CrsMatrix: A class for constructing and using real-valued double-precision sparse compressed r...
const Epetra_Map & DomainMap() const
Returns the Epetra_Map object associated with the domain of this matrix operator. ...
Epetra_CrsMatrix & operator=(const Epetra_CrsMatrix &src)
Assignment operator.
int InputNonlocalValue(int rowoffset, int_type col, double value, int mode)
Epetra_CombineMode
virtual int ReportError(const std::string Message, int ErrorCode) const
Error reporting method.
int N() const
Returns column dimension of system.
const Epetra_CrsGraph & Graph() const
Returns a reference to the Epetra_CrsGraph object associated with this matrix.
Epetra_FECrsMatrix(Epetra_DataAccess CV, const Epetra_Map &RowMap, int *NumEntriesPerRow, bool ignoreNonLocalEntries=false)
Constructor.
Epetra_DataAccess
int * Values()
Returns pointer to the values in vector.
const Epetra_BlockMap & Map() const
Map() method inherited from Epetra_DistObject.
std::vector< double > workData_
int ReplaceGlobalValues(int GlobalRow, int NumEntries, const double *Values, const int *Indices)
override base-class Epetra_CrsMatrix::ReplaceGlobalValues method
virtual ~Epetra_FECrsMatrix()
Destructor.
const Epetra_Comm & Comm() const
Returns a pointer to the Epetra_Comm communicator associated with this matrix.
virtual int ReplaceGlobalValues(int GlobalRow, int NumEntries, const double *Values, const int *Indices)
Replace specified existing values with this list of entries for a given global row of the matrix...
int M() const
Returns row dimension of system.
Epetra_CrsMatrix(Epetra_DataAccess CV, const Epetra_Map &RowMap, const int *NumEntriesPerRow, bool StaticProfile=false)
Epetra_CrsMatrix constructor with variable number of indices per row.
Epetra_CrsGraph: A class for constructing and using sparse compressed row graphs. ...
Epetra_CrsMatrix * nonlocalMatrix_
bool GlobalIndicesTypeMatch(const Epetra_BlockMap &other) const