Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
gpujacobi.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2013-2015: G-CSC, Goethe University Frankfurt
3 * Author: Martin Rupp
4 *
5 * This file is part of UG4.
6 *
7 * UG4 is free software: you can redistribute it and/or modify it under the
8 * terms of the GNU Lesser General Public License version 3 (as published by the
9 * Free Software Foundation) with the following additional attribution
10 * requirements (according to LGPL/GPL v3 §7):
11 *
12 * (1) The following notice must be displayed in the Appropriate Legal Notices
13 * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14 *
15 * (2) The following notice must be displayed at a prominent place in the
16 * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17 *
18 * (3) The following bibliography is recommended for citation and must be
19 * preserved in all covered files:
20 * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21 * parallel geometric multigrid solver on hierarchically distributed grids.
22 * Computing and visualization in science 16, 4 (2013), 151-164"
23 * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24 * flexible software system for simulating pde based models on high performance
25 * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26 *
27 * This program is distributed in the hope that it will be useful,
28 * but WITHOUT ANY WARRANTY; without even the implied warranty of
29 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 * GNU Lesser General Public License for more details.
31 */
32
33#ifdef CUDA_AVAILABLE
34#ifdef UG_GPU
35#ifndef __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__GPUJACOBI__
36#define __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__GPUJACOBI__
37
41
42#ifdef UG_PARALLEL
44#endif
45
49
50
51extern "C" bool
52CUDA_JacobiApply(const double *diagInv, double *corr, const double *defect, const int N);
53
54namespace ug{
55
56
58class GPUJacobi : public IPreconditioner<GPUAlgebra>
59{
60 public:
62 typedef GPUAlgebra TAlgebra;
63 typedef TAlgebra algebra_type;
64
66 typedef TAlgebra::vector_type vector_type;
67
69 typedef TAlgebra::matrix_type matrix_type;
70
72 typedef IPreconditioner<TAlgebra>::matrix_operator_type matrix_operator_type;
73
75 typedef IPreconditioner<TAlgebra> base_type;
76
77 protected:
78 using base_type::set_debug;
79 using base_type::debug_writer;
80 using base_type::write_debug;
81 using base_type::damping;
82
83 public:
85 GPUJacobi() {this->set_damp(1.0);};
86
88 GPUJacobi(number damp) {this->set_damp(damp);};
89
91 GPUJacobi( const GPUJacobi &parent )
92 : base_type(parent)
93 {
94 }
95
98 {
99 return make_sp(new GPUJacobi(*this));
100 }
101
102
104 virtual ~GPUJacobi()
105 {};
106
107 virtual bool supports_parallel() const { return false; }
108 protected:
110 virtual const char* name() const {return "GPUJacobi";}
111
113 virtual bool preprocess(SmartPtr<MatrixOperator<matrix_type, vector_type> > pOp)
114 {
115 PROFILE_BEGIN_GROUP(GPUJacobi_preprocess, "algebra GPUJacobi");
116
117 matrix_type &mat = *pOp;
118 // create help vector to apply diagonal
119 size_t size = mat.num_rows();
120 if(size != mat.num_cols())
121 {
122 UG_LOG("Square Matrix needed for GPUJacobi Iteration.\n");
123 return false;
124 }
125
126 // get damping in constant case to damp at once
127 number damp = 1.0;
128 if(damping()->constant_damping())
129 damp = damping()->damping();
130
131 if(CheckDiagonalInvertible(mat)==false)
132 return false;
133// UG_ASSERT(CheckDiagonalInvertible(mat), "GPUJacobi: A has noninvertible diagonal");
134
135
136 m_diagInv.resize(mat.num_rows());
137 // invert diagonal and multiply by damping
138 for(size_t i = 0; i < mat.num_rows(); ++i)
139 m_diagInv[i] = damp/mat(i,i);
140
141 CUDAHelper::get_instance();
142 // done
143 return true;
144 }
145
146 virtual bool step(SmartPtr<MatrixOperator<matrix_type, vector_type> > pOp, vector_type& c, const vector_type& d)
147 {
148#ifdef UG_PARALLEL
149 // the computed correction is additive
150 c.set_storage_type(PST_ADDITIVE);
151
152 // we make it consistent
153 if(!c.change_storage_type(PST_CONSISTENT))
154 {
155 UG_LOG("ERROR in 'JacobiPreconditioner::apply': "
156 "Cannot change parallel status of correction to consistent.\n");
157 return false;
158 }
159#endif
160
161 const double *devDiagInv = m_diagInv.get_dev_ptr();
162 double *devC = c.get_dev_ptr();
163 const double *devD = d.get_dev_ptr();
164 int N = c.size();
165// UG_LOG("devDiagInv = " << devDiagInv<< " " << m_diagInv.size() << " devC = " << devC << " " << c.size() << " devD = "
166// << devD << " " << d.size() << " N = " << N << "\n");
167 CUDA_JacobiApply(devDiagInv, devC, devD, N);
168
169
170 return true;
171 }
172
174 virtual bool postprocess() {return true;}
175
176 protected:
177 GPUVector<double> m_diagInv;
178
179};
180
181template<>
182class Jacobi<GPUAlgebra> : public GPUJacobi
183{
184 using GPUJacobi::set_damp;
185public:
187 Jacobi() {this->set_damp(1.0);};
188
190 Jacobi(number damp) {this->set_damp(damp);};
191
192 void set_block(bool b)
193 {
194 // always scalar
195 }
196};
197
198
199} // end namespace ug
200
201#endif
202#endif
203#endif // USE_CUDA
location name
Definition checkpoint_util.lua:128
Definition smart_pointer.h:108
void set_damp(SmartPtr< IDamping< TAlgebra::vector_type, TAlgebra::vector_type > > spScaling)
sets a scaling for the correction
Definition linear_iterator.h:163
Jacobi()
default constructor
Definition jacobi.h:117
void set_block(bool b)
sets if blocked jacobi is used (inverting block-diagonal), or plain (scalar) diagonal if false
Definition jacobi.h:145
bool CheckDiagonalInvertible(const TSparseMatrix &A)
Definition sparsematrix_util.h:1094
#define UG_LOG(msg)
Definition log.h:367
double number
Definition types.h:124
CPUAlgebra::matrix_type matrix_type
CPUAlgebra::vector_type vector_type
the ug namespace
#define PROFILE_BEGIN_GROUP(name, groups)
Definition profiler.h:255
SmartPtr< T, FreePolicy > make_sp(T *inst)
returns a SmartPtr for the passed raw pointer
Definition smart_pointer.h:836