pdelabboxassembler.hh 7.92 KB
Newer Older
Andreas Lauser's avatar
Andreas Lauser committed
1
// $Id$
2
3
4
5
6
7
8
9
10
11
12
13
14
15
/*****************************************************************************
 *   Copyright (C) 2009-2010 by Bernd Flemisch                               *
 *   Institute of Hydraulic Engineering                                      *
 *   University of Stuttgart, Germany                                        *
 *   email: <givenname>.<name>@iws.uni-stuttgart.de                          *
 *                                                                           *
 *   This program is free software; you can redistribute it and/or modify    *
 *   it under the terms of the GNU General Public License as published by    *
 *   the Free Software Foundation; either version 2 of the License, or       *
 *   (at your option) any later version, as long as this copyright notice    *
 *   is included in its original form.                                       *
 *                                                                           *
 *   This program is distributed WITHOUT ANY WARRANTY.                       *
 *****************************************************************************/
Andreas Lauser's avatar
Andreas Lauser committed
16
17
#ifndef DUMUX_PDELAB_BOX_ASSEMBLER_HH
#define DUMUX_PDELAB_BOX_ASSEMBLER_HH
18
19
20
21
22
23
24
25

#include<dune/pdelab/finiteelementmap/p1fem.hh>
#include<dune/pdelab/finiteelementmap/q1fem.hh>
#include<dune/pdelab/gridfunctionspace/gridfunctionspace.hh>
#include<dune/pdelab/gridfunctionspace/genericdatahandle.hh>
#include<dune/pdelab/backend/istlvectorbackend.hh>
#include<dune/pdelab/backend/istlmatrixbackend.hh>

Andreas Lauser's avatar
Andreas Lauser committed
26
27
28
29
//#include "pdelabboundarytypes.hh"
#include "pdelabboxlocaloperator.hh"

namespace Dumux {
30

Andreas Lauser's avatar
Andreas Lauser committed
31
namespace PDELab {
32

33
template<class TypeTag>
Andreas Lauser's avatar
Andreas Lauser committed
34
class BoxAssembler
35
{
36
37
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(Model)) Model;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(Problem)) Problem;
38
    enum{numEq = GET_PROP_VALUE(TypeTag, PTAG(NumEq))};
39
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(GridView)) GridView;
40
    enum{dim = GridView::dimension};
41
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(Scalar)) Scalar;
42
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(LocalFEMSpace)) FEM;
Andreas Lauser's avatar
Andreas Lauser committed
43
44
45
46
47
48
49
50

    typedef typename GET_PROP_TYPE(TypeTag, PTAG(Constraints)) Constraints;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(ScalarGridFunctionSpace)) ScalarGridFunctionSpace;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(GridFunctionSpace)) GridFunctionSpace;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(ConstraintsTrafo)) ConstraintsTrafo;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(LocalOperator)) LocalOperator;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(GridOperatorSpace)) GridOperatorSpace;

51
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(LocalJacobian)) LocalJacobian;
Andreas Lauser's avatar
Andreas Lauser committed
52
53
54
55
56
57
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(SolutionVector)) SolutionVector;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(JacobianMatrix)) JacobianMatrix;
public:

    typedef SolutionVector Vector;
    typedef JacobianMatrix Matrix;
58
59
    typedef Matrix RepresentationType;

Andreas Lauser's avatar
Andreas Lauser committed
60
    BoxAssembler()
61
    {
Andreas Lauser's avatar
Andreas Lauser committed
62
        problem_ = 0;
63
64
65
66
67
68
69
70
71
72
        fem_ = 0;
        cn_ = 0;
        scalarGridFunctionSpace_ = 0;
        gridFunctionSpace_ = 0;
        constraintsTrafo_ = 0;
        localOperator_ = 0;
        gridOperatorSpace_ = 0;
        matrix_ = 0;
    }

Andreas Lauser's avatar
Andreas Lauser committed
73
    ~BoxAssembler()
74
75
76
77
78
79
80
81
82
83
84
    {
        delete matrix_;
        delete gridOperatorSpace_;
        delete localOperator_;
        delete constraintsTrafo_;
        delete gridFunctionSpace_;
        delete scalarGridFunctionSpace_;
        delete cn_;
        delete fem_;
    }

Andreas Lauser's avatar
Andreas Lauser committed
85
86
87
88
89
90
91
92
93
94
    void init(Problem& problem)
    {
        problem_ = &problem;
        fem_ = new FEM();
        //cn_ = new Constraints(*problem_);
        cn_ = new Constraints();
        scalarGridFunctionSpace_ = new ScalarGridFunctionSpace(problem_->gridView(), *fem_, *cn_);
        gridFunctionSpace_ = new GridFunctionSpace(*scalarGridFunctionSpace_);

        //cn_->compute_ghosts(*gridFunctionSpace_);
95

Andreas Lauser's avatar
Andreas Lauser committed
96
97
98
99
100
101
        //typedef BoundaryIndexHelper<TypeTag> BoundaryFunction;
        //BoundaryFunction *bTypes = new BoundaryFunction();
        constraintsTrafo_ = new ConstraintsTrafo();
        //Dune::PDELab::constraints(*bTypes, *gridFunctionSpace_, *constraintsTrafo_, false);

        localOperator_ = new LocalOperator(problem_->model());
102
        gridOperatorSpace_ =
Andreas Lauser's avatar
Andreas Lauser committed
103
104
105
106
107
            new GridOperatorSpace(*gridFunctionSpace_, *constraintsTrafo_,
                                  *gridFunctionSpace_, *constraintsTrafo_, *localOperator_);

        matrix_ = new Matrix(*gridOperatorSpace_);
        *matrix_ = 0;
108
        reuseMatrix_ = false;
Andreas Lauser's avatar
Andreas Lauser committed
109
110
    }

111
112
    void assemble(SolutionVector &u)
    {
113
114
115
116
117
        if (!reuseMatrix_) {
            *matrix_ = 0;
            gridOperatorSpace_->jacobian(u, *matrix_);
        }
        reuseMatrix_ = false;
118
119
120
121
122
123
124
125
126
127

        residual_.base().resize(u.size());
        residual_ = 0;
        gridOperatorSpace_->residual(u, residual_);

#if 0
        // rescale jacobian and right hand side to the largest
        // entry on the main diagonal block matrix
        typedef typename Matrix::RowIterator RowIterator;
        typedef typename Matrix::ColIterator ColIterator;
128
        typedef typename Matrix::block_type BlockType;
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
        const typename Matrix::block_type::size_type rowsInBlock = Matrix::block_type::rows;
        const typename Matrix::block_type::size_type colsInBlock = Matrix::block_type::cols;
        Scalar diagonalEntry[rowsInBlock];
        Vector diagonalEntries(*f);
        RowIterator endIBlock = matrix_->end();
        for (RowIterator iBlock = matrix_->begin(); iBlock != endIBlock; ++iBlock) {
            BlockType &diagBlock = (*iBlock)[iBlock.index()];

            for (int i = 0; i < rowsInBlock; ++i) {
                diagonalEntry[i] = 0;
                for (int j = 0; j < colsInBlock; ++j) {
                    diagonalEntry[i] = std::max(diagonalEntry[i],
                                                std::abs(diagBlock[i][j]));
                }

                if (diagonalEntry[i] < 1e-14)
                    diagonalEntry[i] = 1.0;

                diagonalEntries[iBlock.index()][i] = diagonalEntry[i];
            }
        }

        Dune::PDELab::AddDataHandle<GridFunctionSpace,Vector> adddh(*gridFunctionSpace_, diagonalEntries);
        if (gridFunctionSpace_->gridview().comm().size()>1)
          gridFunctionSpace_->gridview().communicate(adddh,
              Dune::InteriorBorder_InteriorBorder_Interface,
              Dune::ForwardCommunication);

        for (RowIterator iBlock = matrix_->begin(); iBlock != endIBlock; ++iBlock) {
            // divide right-hand side
            for (int i = 0; i < rowsInBlock; i++) {
                (*f)[iBlock.index()][i] /= diagonalEntries[iBlock.index()][i];
            }

            // divide row of the jacobian
            ColIterator endJBlock = iBlock->end();
            for (ColIterator jBlock = iBlock->begin(); jBlock != endJBlock; ++jBlock) {
                for (int i = 0; i < rowsInBlock; i++) {
                    for (int j = 0; j < colsInBlock; j++) {
                        (*jBlock)[i][j] /= diagonalEntries[iBlock.index()][i];
                    }
                }
            }
        }
#endif
    }

176
177
178
179
    void setMatrixReuseable(bool yesno = true)
    { reuseMatrix_ = yesno; }
    

180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
    const GridFunctionSpace& gridFunctionSpace() const
    {
        return *gridFunctionSpace_;
    }

    const ConstraintsTrafo& constraintsTrafo() const
    {
        return *constraintsTrafo_;
    }

    //! return const reference to matrix
    const Matrix& matrix() const
    { return *matrix_; }

    //! return const reference to residual
    const SolutionVector& residual() const
    { return residual_; }

private:
Andreas Lauser's avatar
Andreas Lauser committed
199
    Problem *problem_;
200
201
202
203
204
205
206
207
208
    Constraints *cn_;
    FEM *fem_;
    ScalarGridFunctionSpace *scalarGridFunctionSpace_;
    GridFunctionSpace *gridFunctionSpace_;
    ConstraintsTrafo *constraintsTrafo_;
    LocalOperator *localOperator_;
    GridOperatorSpace *gridOperatorSpace_;

    Matrix *matrix_;
209
    bool reuseMatrix_;
210
211
212
    SolutionVector residual_;
};

Andreas Lauser's avatar
Andreas Lauser committed
213
214
} // namespace PDELab
} // namespace Dumux
215
216

#endif