EpetraExt
Development
Toggle main menu visibility
Loading...
Searching...
No Matches
src
transform
EpetraExt_Reindex_LinearProblem.cpp
Go to the documentation of this file.
1
//@HEADER
2
// ***********************************************************************
3
//
4
// EpetraExt: Epetra Extended - Linear Algebra Services Package
5
// Copyright (2011) Sandia Corporation
6
//
7
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8
// the U.S. Government retains certain rights in this software.
9
//
10
// Redistribution and use in source and binary forms, with or without
11
// modification, are permitted provided that the following conditions are
12
// met:
13
//
14
// 1. Redistributions of source code must retain the above copyright
15
// notice, this list of conditions and the following disclaimer.
16
//
17
// 2. Redistributions in binary form must reproduce the above copyright
18
// notice, this list of conditions and the following disclaimer in the
19
// documentation and/or other materials provided with the distribution.
20
//
21
// 3. Neither the name of the Corporation nor the names of the
22
// contributors may be used to endorse or promote products derived from
23
// this software without specific prior written permission.
24
//
25
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36
//
37
// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38
//
39
// ***********************************************************************
40
//@HEADER
41
42
#include <
EpetraExt_Reindex_LinearProblem.h
>
43
44
#include <
EpetraExt_Reindex_CrsMatrix.h
>
45
#include <
EpetraExt_Reindex_MultiVector.h
>
46
47
#include <Epetra_LinearProblem.h>
48
#include <Epetra_CrsMatrix.h>
49
#include <Epetra_MultiVector.h>
50
#include <Epetra_Map.h>
51
52
namespace
EpetraExt
{
53
54
LinearProblem_Reindex::
55
~LinearProblem_Reindex
()
56
{
57
if
(
newObj_
)
delete
newObj_
;
58
59
if
( MatTrans_ )
delete
MatTrans_;
60
if
( LHSTrans_ )
delete
LHSTrans_;
61
if
( RHSTrans_ )
delete
RHSTrans_;
62
63
if
( NewRowMapOwned_ )
delete
NewRowMap_;
64
}
65
66
LinearProblem_Reindex::NewTypeRef
67
LinearProblem_Reindex::
68
operator()
(
OriginalTypeRef
orig )
69
{
70
Epetra_CrsMatrix
* OldMatrix =
dynamic_cast<
Epetra_CrsMatrix
*
>
( orig.GetMatrix() );
71
Epetra_MultiVector
* OldRHS = orig.GetRHS();
72
Epetra_MultiVector
* OldLHS = orig.GetLHS();
73
const
Epetra_BlockMap
& OldRowMap = OldMatrix->
Map
();
74
75
//If new map not passed in create one with lex ordering
76
if
( !NewRowMap_ )
77
{
78
int
NumMyElements = OldRowMap.
NumMyElements
();
79
80
#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
81
if
(OldRowMap.
GlobalIndicesInt
()) {
82
int
NumGlobalElements = OldRowMap.
NumGlobalElements
();
83
NewRowMap_ =
new
Epetra_Map
( NumGlobalElements, NumMyElements, 0, OldRowMap.
Comm
() );
84
}
85
else
86
#endif
87
#ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
88
if
(OldRowMap.
GlobalIndicesLongLong
()) {
89
long
long
NumGlobalElements = OldRowMap.NumGlobalElements64();
90
NewRowMap_ =
new
Epetra_Map
( NumGlobalElements, NumMyElements, 0LL, OldRowMap.
Comm
() );
91
}
92
else
93
#endif
94
throw
"LinearProblem_Reindex::operator(): GlobalIndices type unknown for OldRowMap"
;
95
96
NewRowMapOwned_ =
true
;
97
}
98
99
MatTrans_ =
new
CrsMatrix_Reindex
( *NewRowMap_ );
100
LHSTrans_ =
new
MultiVector_Reindex
( *NewRowMap_ );
101
RHSTrans_ =
new
MultiVector_Reindex
( *NewRowMap_ );
102
103
Epetra_CrsMatrix
* NewMatrix = &((*MatTrans_)( *OldMatrix ));
104
Epetra_MultiVector
* NewLHS = &((*LHSTrans_)( *OldLHS ));
105
Epetra_MultiVector
* NewRHS = &((*RHSTrans_)( *OldRHS ));
106
107
newObj_
=
new
Epetra_LinearProblem
( NewMatrix, NewLHS, NewRHS );
108
109
return
*
newObj_
;
110
}
111
112
}
//namespace EpetraExt
113
EpetraExt_Reindex_CrsMatrix.h
EpetraExt_Reindex_LinearProblem.h
EpetraExt_Reindex_MultiVector.h
EpetraExt::CrsMatrix_Reindex
Given an Epetra_CrsMatrix, a "reindexed" version is returned based on the new row map.
Definition
EpetraExt_Reindex_CrsMatrix.h:57
EpetraExt::LinearProblem_Reindex::operator()
NewTypeRef operator()(OriginalTypeRef orig)
Constructs a new view the original LP, "reindexed" using the given NewRowMap.
Definition
EpetraExt_Reindex_LinearProblem.cpp:68
EpetraExt::LinearProblem_Reindex::~LinearProblem_Reindex
~LinearProblem_Reindex()
Destructor.
Definition
EpetraExt_Reindex_LinearProblem.cpp:55
EpetraExt::MultiVector_Reindex
Given an input Epetra_MultiVector, a "reindexed" view is returned.
Definition
EpetraExt_Reindex_MultiVector.h:55
EpetraExt::Transform< Epetra_LinearProblem, Epetra_LinearProblem >::OriginalTypeRef
Epetra_LinearProblem & OriginalTypeRef
Definition
EpetraExt_Transform.h:74
EpetraExt::Transform< Epetra_LinearProblem, Epetra_LinearProblem >::newObj_
NewTypePtr newObj_
Definition
EpetraExt_Transform.h:218
EpetraExt::Transform< Epetra_LinearProblem, Epetra_LinearProblem >::NewTypeRef
Epetra_LinearProblem & NewTypeRef
Definition
EpetraExt_Transform.h:79
Epetra_BlockMap
Epetra_BlockMap::GlobalIndicesInt
bool GlobalIndicesInt() const
Epetra_BlockMap::NumGlobalElements
int NumGlobalElements() const
Epetra_BlockMap::Comm
const Epetra_Comm & Comm() const
Epetra_BlockMap::NumMyElements
int NumMyElements() const
Epetra_BlockMap::GlobalIndicesLongLong
bool GlobalIndicesLongLong() const
Epetra_CrsMatrix
Epetra_CrsMatrix::Map
const Epetra_BlockMap & Map() const
Epetra_LinearProblem
Epetra_Map
Epetra_MultiVector
EpetraExt
EpetraExt::BlockCrsMatrix: A class for constructing a distributed block matrix.
Definition
EpetraExt_BlockCrsMatrix.cpp:46
Generated by
1.17.0