-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathbatched_8hpp_source.html
317 lines (315 loc) · 53.6 KB
/
batched_8hpp_source.html
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
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
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
<!-- HTML header for doxygen 1.9.6-->
<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "https://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
<html xmlns="http://www.w3.org/1999/xhtml" lang="en-US">
<head>
<meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/>
<meta http-equiv="X-UA-Compatible" content="IE=11"/>
<meta name="generator" content="Doxygen 1.11.0"/>
<meta name="viewport" content="width=device-width, initial-scale=1"/>
<title>MFEM: linalg/batched/batched.hpp Source File</title>
<link href="tabs.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="jquery.js"></script>
<script type="text/javascript" src="dynsections.js"></script>
<link href="navtree.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="resize.js"></script>
<script type="text/javascript" src="cookie.js"></script>
<link href="search/search.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="search/searchdata.js"></script>
<script type="text/javascript" src="search/search.js"></script>
<script type="text/javascript">
window.MathJax = {
options: {
ignoreHtmlClass: 'tex2jax_ignore',
processHtmlClass: 'tex2jax_process'
},
loader: {
load: ['[tex]/ams']
},
tex: {
macros: {},
packages: ['base','configmacros','ams']
}
};
</script>
<script type="text/javascript" id="MathJax-script" async="async" src="https://cdn.jsdelivr.net/npm/mathjax@3/es5/tex-chtml.js"></script>
<link href="doxygen.css" rel="stylesheet" type="text/css" />
<link href="doxygen-awesome.css" rel="stylesheet" type="text/css"/>
<link href="customization.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="doxygen-awesome-darkmode-toggle.js"></script>
<script type="text/javascript">
DoxygenAwesomeDarkModeToggle.init()
</script>
</head>
<body>
<div id="top"><!-- do not remove this div, it is closed by doxygen! -->
<div id="titlearea">
<table cellspacing="0" cellpadding="0">
<tbody>
<tr id="projectrow">
<td id="projectlogo"><img alt="Logo" src="logo-small.png"/></td>
<td id="projectalign">
<div id="projectname">MFEM<span id="projectnumber"> v4.8.0</span>
</div>
<div id="projectbrief">Finite element discretization library</div>
</td>
</tr>
</tbody>
</table>
</div>
<!-- end header part -->
<!-- Generated by Doxygen 1.11.0 -->
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699&dn=expat.txt MIT */
var searchBox = new SearchBox("searchBox", "search/",'.html');
/* @license-end */
</script>
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699&dn=expat.txt MIT */
$(function() { codefold.init(0); });
/* @license-end */
</script>
<script type="text/javascript" src="menudata.js"></script>
<script type="text/javascript" src="menu.js"></script>
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699&dn=expat.txt MIT */
$(function() {
initMenu('',true,false,'search.php','Search',false);
$(function() { init_search(); });
});
/* @license-end */
</script>
<div id="main-nav"></div>
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699&dn=expat.txt MIT */
$(function(){ initResizable(false); });
/* @license-end */
</script>
<!-- window showing the filter options -->
<div id="MSearchSelectWindow"
onmouseover="return searchBox.OnSearchSelectShow()"
onmouseout="return searchBox.OnSearchSelectHide()"
onkeydown="return searchBox.OnSearchSelectKey(event)">
</div>
<!-- iframe showing the search results (closed by default) -->
<div id="MSearchResultsWindow">
<div id="MSearchResults">
<div class="SRPage">
<div id="SRIndex">
<div id="SRResults"></div>
<div class="SRStatus" id="Loading">Loading...</div>
<div class="SRStatus" id="Searching">Searching...</div>
<div class="SRStatus" id="NoMatches">No Matches</div>
</div>
</div>
</div>
</div>
<div id="nav-path" class="navpath">
<ul>
<li class="navelem"><a class="el" href="dir_c1845e5bf011c2136883d0d7188f50e9.html">linalg</a></li><li class="navelem"><a class="el" href="dir_7c63c39fa17c0116e286dd47ad28981a.html">batched</a></li> </ul>
</div>
</div><!-- top -->
<div id="doc-content">
<div class="header">
<div class="headertitle"><div class="title">batched.hpp</div></div>
</div><!--header-->
<div class="contents">
<a href="batched_8hpp.html">Go to the documentation of this file.</a><div class="fragment"><div class="line"><a id="l00001" name="l00001"></a><span class="lineno"> 1</span><span class="comment">// Copyright (c) 2010-2025, Lawrence Livermore National Security, LLC. Produced</span></div>
<div class="line"><a id="l00002" name="l00002"></a><span class="lineno"> 2</span><span class="comment">// at the Lawrence Livermore National Laboratory. All Rights reserved. See files</span></div>
<div class="line"><a id="l00003" name="l00003"></a><span class="lineno"> 3</span><span class="comment">// LICENSE and NOTICE for details. LLNL-CODE-806117.</span></div>
<div class="line"><a id="l00004" name="l00004"></a><span class="lineno"> 4</span><span class="comment">//</span></div>
<div class="line"><a id="l00005" name="l00005"></a><span class="lineno"> 5</span><span class="comment">// This file is part of the MFEM library. For more information and source code</span></div>
<div class="line"><a id="l00006" name="l00006"></a><span class="lineno"> 6</span><span class="comment">// availability visit https://mfem.org.</span></div>
<div class="line"><a id="l00007" name="l00007"></a><span class="lineno"> 7</span><span class="comment">//</span></div>
<div class="line"><a id="l00008" name="l00008"></a><span class="lineno"> 8</span><span class="comment">// MFEM is free software; you can redistribute it and/or modify it under the</span></div>
<div class="line"><a id="l00009" name="l00009"></a><span class="lineno"> 9</span><span class="comment">// terms of the BSD-3 license. We welcome feedback and contributions, see file</span></div>
<div class="line"><a id="l00010" name="l00010"></a><span class="lineno"> 10</span><span class="comment">// CONTRIBUTING.md for details.</span></div>
<div class="line"><a id="l00011" name="l00011"></a><span class="lineno"> 11</span> </div>
<div class="line"><a id="l00012" name="l00012"></a><span class="lineno"> 12</span><span class="preprocessor">#ifndef MFEM_BATCHED_LINALG</span></div>
<div class="line"><a id="l00013" name="l00013"></a><span class="lineno"> 13</span><span class="preprocessor">#define MFEM_BATCHED_LINALG</span></div>
<div class="line"><a id="l00014" name="l00014"></a><span class="lineno"> 14</span> </div>
<div class="line"><a id="l00015" name="l00015"></a><span class="lineno"> 15</span><span class="preprocessor">#include "<a class="code" href="config_8hpp.html">../../config/config.hpp</a>"</span></div>
<div class="line"><a id="l00016" name="l00016"></a><span class="lineno"> 16</span><span class="preprocessor">#include "<a class="code" href="densemat_8hpp.html">../densemat.hpp</a>"</span></div>
<div class="line"><a id="l00017" name="l00017"></a><span class="lineno"> 17</span><span class="preprocessor">#include <array></span></div>
<div class="line"><a id="l00018" name="l00018"></a><span class="lineno"> 18</span><span class="preprocessor">#include <memory></span></div>
<div class="line"><a id="l00019" name="l00019"></a><span class="lineno"> 19</span> </div>
<div class="line"><a id="l00020" name="l00020"></a><span class="lineno"> 20</span><span class="keyword">namespace </span><a class="code hl_namespace" href="namespacemfem.html">mfem</a></div>
<div class="line"><a id="l00021" name="l00021"></a><span class="lineno"> 21</span>{</div>
<div class="line"><a id="l00022" name="l00022"></a><span class="lineno"> 22</span><span class="comment"></span> </div>
<div class="line"><a id="l00023" name="l00023"></a><span class="lineno"> 23</span><span class="comment">/// @brief Class for performing batched linear algebra operations, potentially</span></div>
<div class="line"><a id="l00024" name="l00024"></a><span class="lineno"> 24</span><span class="comment">/// using accelerated algorithms (GPU BLAS or MAGMA). Accessed using static</span></div>
<div class="line"><a id="l00025" name="l00025"></a><span class="lineno"> 25</span><span class="comment">/// member functions.</span></div>
<div class="line"><a id="l00026" name="l00026"></a><span class="lineno"> 26</span><span class="comment">///</span></div>
<div class="line"><a id="l00027" name="l00027"></a><span class="lineno"> 27</span><span class="comment">/// The static member functions will delegate to the active backend (which can</span></div>
<div class="line"><a id="l00028" name="l00028"></a><span class="lineno"> 28</span><span class="comment">/// be set using SetActiveBackend(), see BatchedLinAlg::Backend for all</span></div>
<div class="line"><a id="l00029" name="l00029"></a><span class="lineno"> 29</span><span class="comment">/// available backends and the order in which they will be chosen initially).</span></div>
<div class="line"><a id="l00030" name="l00030"></a><span class="lineno"> 30</span><span class="comment">/// Operations can be performed directly with a specific backend using Get().</span></div>
<div class="foldopen" id="foldopen00031" data-start="{" data-end="};">
<div class="line"><a id="l00031" name="l00031"></a><span class="lineno"><a class="line" href="classmfem_1_1BatchedLinAlg.html"> 31</a></span><span class="comment"></span><span class="keyword">class </span><a class="code hl_class" href="classmfem_1_1BatchedLinAlg.html">BatchedLinAlg</a></div>
<div class="line"><a id="l00032" name="l00032"></a><span class="lineno"> 32</span>{</div>
<div class="line"><a id="l00033" name="l00033"></a><span class="lineno"> 33</span><span class="keyword">public</span>:<span class="comment"></span></div>
<div class="line"><a id="l00034" name="l00034"></a><span class="lineno"> 34</span><span class="comment"> /// @brief Available backends for implementations of batched algorithms.</span></div>
<div class="line"><a id="l00035" name="l00035"></a><span class="lineno"> 35</span><span class="comment"> ///</span></div>
<div class="line"><a id="l00036" name="l00036"></a><span class="lineno"> 36</span><span class="comment"> /// The initially active backend will be the first available backend in this</span></div>
<div class="line"><a id="l00037" name="l00037"></a><span class="lineno"> 37</span><span class="comment"> /// order: MAGMA, GPU_BLAS, NATIVE.</span></div>
<div class="foldopen" id="foldopen00038" data-start="{" data-end="};">
<div class="line"><a id="l00038" name="l00038"></a><span class="lineno"><a class="line" href="classmfem_1_1BatchedLinAlg.html#a40da089d2e89a7c7ecef22423473b6ca"> 38</a></span><span class="comment"></span> <span class="keyword">enum</span> <a class="code hl_enumeration" href="classmfem_1_1BatchedLinAlg.html#a40da089d2e89a7c7ecef22423473b6ca">Backend</a></div>
<div class="line"><a id="l00039" name="l00039"></a><span class="lineno"> 39</span> {<span class="comment"></span></div>
<div class="line"><a id="l00040" name="l00040"></a><span class="lineno"> 40</span><span class="comment"> /// @brief The standard MFEM backend, implemented using mfem::forall</span></div>
<div class="line"><a id="l00041" name="l00041"></a><span class="lineno"> 41</span><span class="comment"> /// kernels. Not as performant as the other kernels.</span></div>
<div class="line"><a id="l00042" name="l00042"></a><span class="lineno"><a class="line" href="classmfem_1_1BatchedLinAlg.html#a40da089d2e89a7c7ecef22423473b6caadd4da6e88f75868880f6692ed6b5eeb4"> 42</a></span><span class="comment"></span> <a class="code hl_enumvalue" href="classmfem_1_1BatchedLinAlg.html#a40da089d2e89a7c7ecef22423473b6caadd4da6e88f75868880f6692ed6b5eeb4">NATIVE</a>,<span class="comment"></span></div>
<div class="line"><a id="l00043" name="l00043"></a><span class="lineno"> 43</span><span class="comment"> /// @brief Either cuBLAS or hipBLAS, depending on whether MFEM is using</span></div>
<div class="line"><a id="l00044" name="l00044"></a><span class="lineno"> 44</span><span class="comment"> /// CUDA or HIP. Not available otherwise.</span></div>
<div class="line"><a id="l00045" name="l00045"></a><span class="lineno"><a class="line" href="classmfem_1_1BatchedLinAlg.html#a40da089d2e89a7c7ecef22423473b6caa388e52e74ac0158131847fcc330b16ee"> 45</a></span><span class="comment"></span> <a class="code hl_enumvalue" href="classmfem_1_1BatchedLinAlg.html#a40da089d2e89a7c7ecef22423473b6caa388e52e74ac0158131847fcc330b16ee">GPU_BLAS</a>,<span class="comment"></span></div>
<div class="line"><a id="l00046" name="l00046"></a><span class="lineno"> 46</span><span class="comment"> /// MAGMA backend, only available if MFEM is compiled with MAGMA support.</span></div>
<div class="line"><a id="l00047" name="l00047"></a><span class="lineno"><a class="line" href="classmfem_1_1BatchedLinAlg.html#a40da089d2e89a7c7ecef22423473b6caa462cb27cf8966b995bce9663e7d327ce"> 47</a></span><span class="comment"></span> <a class="code hl_enumvalue" href="classmfem_1_1BatchedLinAlg.html#a40da089d2e89a7c7ecef22423473b6caa462cb27cf8966b995bce9663e7d327ce">MAGMA</a>,<span class="comment"></span></div>
<div class="line"><a id="l00048" name="l00048"></a><span class="lineno"> 48</span><span class="comment"> /// Counter for the number of backends.</span></div>
<div class="line"><a id="l00049" name="l00049"></a><span class="lineno"> 49</span><span class="comment"></span> <a class="code hl_enumvalue" href="classmfem_1_1BatchedLinAlg.html#a40da089d2e89a7c7ecef22423473b6caa467cc6ba51d7366f915ac4413ab5ef09">NUM_BACKENDS</a></div>
<div class="line"><a id="l00050" name="l00050"></a><span class="lineno"><a class="line" href="classmfem_1_1BatchedLinAlg.html#a40da089d2e89a7c7ecef22423473b6caa467cc6ba51d7366f915ac4413ab5ef09"> 50</a></span> };</div>
</div>
<div class="line"><a id="l00051" name="l00051"></a><span class="lineno"> 51</span><span class="comment"></span> </div>
<div class="line"><a id="l00052" name="l00052"></a><span class="lineno"> 52</span><span class="comment"> /// Operation type (transposed or not transposed)</span></div>
<div class="foldopen" id="foldopen00053" data-start="{" data-end="};">
<div class="line"><a id="l00053" name="l00053"></a><span class="lineno"><a class="line" href="classmfem_1_1BatchedLinAlg.html#ade3f74dd4281b3e06ab1cee7448b0d8d"> 53</a></span><span class="comment"></span> <span class="keyword">enum</span> <a class="code hl_enumeration" href="classmfem_1_1BatchedLinAlg.html#ade3f74dd4281b3e06ab1cee7448b0d8d">Op</a></div>
<div class="line"><a id="l00054" name="l00054"></a><span class="lineno"> 54</span> {</div>
<div class="line"><a id="l00055" name="l00055"></a><span class="lineno"><a class="line" href="classmfem_1_1BatchedLinAlg.html#ade3f74dd4281b3e06ab1cee7448b0d8dabefd6c3489b5cd88b472b82cff88df8f"> 55</a></span> <a class="code hl_enumvalue" href="classmfem_1_1BatchedLinAlg.html#ade3f74dd4281b3e06ab1cee7448b0d8dabefd6c3489b5cd88b472b82cff88df8f">N</a>, <span class="comment">///< Not transposed.</span></div>
<div class="line"><a id="l00056" name="l00056"></a><span class="lineno"> 56</span> <a class="code hl_enumvalue" href="classmfem_1_1BatchedLinAlg.html#ade3f74dd4281b3e06ab1cee7448b0d8da19710f7a1d92fa9659d19cef79fe781b">T</a> <span class="comment">///< Transposed.</span></div>
<div class="line"><a id="l00057" name="l00057"></a><span class="lineno"><a class="line" href="classmfem_1_1BatchedLinAlg.html#ade3f74dd4281b3e06ab1cee7448b0d8da19710f7a1d92fa9659d19cef79fe781b"> 57</a></span> };</div>
</div>
<div class="line"><a id="l00058" name="l00058"></a><span class="lineno"> 58</span> </div>
<div class="line"><a id="l00059" name="l00059"></a><span class="lineno"> 59</span><span class="keyword">private</span>:<span class="comment"></span></div>
<div class="line"><a id="l00060" name="l00060"></a><span class="lineno"> 60</span><span class="comment"> /// All available backends. Unavailable backends will be nullptr.</span></div>
<div class="line"><a id="l00061" name="l00061"></a><span class="lineno"> 61</span><span class="comment"></span> std::array<std::unique_ptr<class BatchedLinAlgBase>,</div>
<div class="line"><a id="l00062" name="l00062"></a><span class="lineno"> 62</span> <a class="code hl_enumvalue" href="classmfem_1_1BatchedLinAlg.html#a40da089d2e89a7c7ecef22423473b6caa467cc6ba51d7366f915ac4413ab5ef09">Backend::NUM_BACKENDS</a>> backends;</div>
<div class="line"><a id="l00063" name="l00063"></a><span class="lineno"> 63</span> <a class="code hl_struct" href="structmfem_1_1Backend.html">Backend</a> active_backend;<span class="comment"></span></div>
<div class="line"><a id="l00064" name="l00064"></a><span class="lineno"> 64</span><span class="comment"> /// Default constructor. Private.</span></div>
<div class="line"><a id="l00065" name="l00065"></a><span class="lineno"> 65</span><span class="comment"></span> <a class="code hl_class" href="classmfem_1_1BatchedLinAlg.html">BatchedLinAlg</a>();<span class="comment"></span></div>
<div class="line"><a id="l00066" name="l00066"></a><span class="lineno"> 66</span><span class="comment"> /// Return the singleton instance.</span></div>
<div class="line"><a id="l00067" name="l00067"></a><span class="lineno"> 67</span><span class="comment"></span> <span class="keyword">static</span> <a class="code hl_class" href="classmfem_1_1BatchedLinAlg.html">BatchedLinAlg</a> &Instance();</div>
<div class="line"><a id="l00068" name="l00068"></a><span class="lineno"> 68</span><span class="keyword">public</span>:<span class="comment"></span></div>
<div class="line"><a id="l00069" name="l00069"></a><span class="lineno"> 69</span><span class="comment"> /// @brief Computes $y = \alpha A^{op} x + \beta y$.</span></div>
<div class="line"><a id="l00070" name="l00070"></a><span class="lineno"> 70</span><span class="comment"> ///</span></div>
<div class="line"><a id="l00071" name="l00071"></a><span class="lineno"> 71</span><span class="comment"> /// $A^{op}$ is either $A$ or $A^T$ depending on the value of @a op.</span></div>
<div class="line"><a id="l00072" name="l00072"></a><span class="lineno"> 72</span><span class="comment"> /// $A$ is a block diagonal matrix, represented by the DenseTensor @a A with</span></div>
<div class="line"><a id="l00073" name="l00073"></a><span class="lineno"> 73</span><span class="comment"> /// shape (m, n, n_mat). $x$ has shape (tr?m:n, k, n_mat), and $y$ has shape</span></div>
<div class="line"><a id="l00074" name="l00074"></a><span class="lineno"> 74</span><span class="comment"> /// (tr?n:m, k, n_mat), where 'tr' is true in the transposed case.</span></div>
<div class="line"><a id="l00075" name="l00075"></a><span class="lineno"> 75</span><span class="comment"></span> <span class="keyword">static</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classmfem_1_1BatchedLinAlg.html#ab3676ca00da6882e00dfb7718b5241dd">AddMult</a>(<span class="keyword">const</span> <a class="code hl_class" href="classmfem_1_1DenseTensor.html">DenseTensor</a> &A, <span class="keyword">const</span> <a class="code hl_class" href="classmfem_1_1Vector.html">Vector</a> &x, <a class="code hl_class" href="classmfem_1_1Vector.html">Vector</a> &y,</div>
<div class="line"><a id="l00076" name="l00076"></a><span class="lineno"> 76</span> <a class="code hl_typedef" href="namespacemfem.html#ac8ba85124b5de90bf368176b69d6019c">real_t</a> <a class="code hl_variable" href="ex15_8cpp.html#a7dc6cd505f60a008e9246be53a4d460d">alpha</a> = 1.0, <a class="code hl_typedef" href="namespacemfem.html#ac8ba85124b5de90bf368176b69d6019c">real_t</a> <a class="code hl_variable" href="convection-diffusion_8cpp.html#a20d35965a489ee0e4396d34b3867aaed">beta</a> = 1.0,</div>
<div class="line"><a id="l00077" name="l00077"></a><span class="lineno"> 77</span> <a class="code hl_enumeration" href="classmfem_1_1BatchedLinAlg.html#ade3f74dd4281b3e06ab1cee7448b0d8d">Op</a> op = <a class="code hl_enumvalue" href="classmfem_1_1BatchedLinAlg.html#ade3f74dd4281b3e06ab1cee7448b0d8dabefd6c3489b5cd88b472b82cff88df8f">Op::N</a>);<span class="comment"></span></div>
<div class="line"><a id="l00078" name="l00078"></a><span class="lineno"> 78</span><span class="comment"> /// Computes $y = A x$ (e.g. by calling @ref AddMult "AddMult(A,x,y,1,0,Op::N)").</span></div>
<div class="line"><a id="l00079" name="l00079"></a><span class="lineno"> 79</span><span class="comment"></span> <span class="keyword">static</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classmfem_1_1BatchedLinAlg.html#a327ee59ae1496af66d18207a94d491fd">Mult</a>(<span class="keyword">const</span> <a class="code hl_class" href="classmfem_1_1DenseTensor.html">DenseTensor</a> &A, <span class="keyword">const</span> <a class="code hl_class" href="classmfem_1_1Vector.html">Vector</a> &x, <a class="code hl_class" href="classmfem_1_1Vector.html">Vector</a> &y);<span class="comment"></span></div>
<div class="line"><a id="l00080" name="l00080"></a><span class="lineno"> 80</span><span class="comment"> /// Computes $y = A^T x$ (e.g. by calling @ref AddMult "AddMult(A,x,y,1,0,Op::T)").</span></div>
<div class="line"><a id="l00081" name="l00081"></a><span class="lineno"> 81</span><span class="comment"></span> <span class="keyword">static</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classmfem_1_1BatchedLinAlg.html#a806592aaaa535f422453f75e855735ea">MultTranspose</a>(<span class="keyword">const</span> <a class="code hl_class" href="classmfem_1_1DenseTensor.html">DenseTensor</a> &A, <span class="keyword">const</span> <a class="code hl_class" href="classmfem_1_1Vector.html">Vector</a> &x, <a class="code hl_class" href="classmfem_1_1Vector.html">Vector</a> &y);<span class="comment"></span></div>
<div class="line"><a id="l00082" name="l00082"></a><span class="lineno"> 82</span><span class="comment"> /// @brief Replaces the block diagonal matrix $A$ with its inverse $A^{-1}$.</span></div>
<div class="line"><a id="l00083" name="l00083"></a><span class="lineno"> 83</span><span class="comment"> ///</span></div>
<div class="line"><a id="l00084" name="l00084"></a><span class="lineno"> 84</span><span class="comment"> /// $A$ is represented by the DenseTensor @a A with shape (m, m, n_mat).</span></div>
<div class="line"><a id="l00085" name="l00085"></a><span class="lineno"> 85</span><span class="comment"></span> <span class="keyword">static</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classmfem_1_1BatchedLinAlg.html#ae64e57da2451c570f6e0fa06ebd8a111">Invert</a>(<a class="code hl_class" href="classmfem_1_1DenseTensor.html">DenseTensor</a> &A);<span class="comment"></span></div>
<div class="line"><a id="l00086" name="l00086"></a><span class="lineno"> 86</span><span class="comment"> /// @brief Replaces the block diagonal matrix $A$ with its LU factors. The</span></div>
<div class="line"><a id="l00087" name="l00087"></a><span class="lineno"> 87</span><span class="comment"> /// pivots are stored in @a P.</span></div>
<div class="line"><a id="l00088" name="l00088"></a><span class="lineno"> 88</span><span class="comment"> ///</span></div>
<div class="line"><a id="l00089" name="l00089"></a><span class="lineno"> 89</span><span class="comment"> /// $A$ is represented by the DenseTensor @a A with shape (n, n, n_mat). On</span></div>
<div class="line"><a id="l00090" name="l00090"></a><span class="lineno"> 90</span><span class="comment"> /// output, $P$ has shape (n, n_mat).</span></div>
<div class="line"><a id="l00091" name="l00091"></a><span class="lineno"> 91</span><span class="comment"></span> <span class="keyword">static</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classmfem_1_1BatchedLinAlg.html#aaf20d258a08692f4e9c7e3125d980658">LUFactor</a>(<a class="code hl_class" href="classmfem_1_1DenseTensor.html">DenseTensor</a> &A, <a class="code hl_class" href="classmfem_1_1Array.html">Array<int></a> &P);<span class="comment"></span></div>
<div class="line"><a id="l00092" name="l00092"></a><span class="lineno"> 92</span><span class="comment"> /// @brief Replaces $x$ with $A^{-1} x$, given the LU factors @a A and pivots</span></div>
<div class="line"><a id="l00093" name="l00093"></a><span class="lineno"> 93</span><span class="comment"> /// @a P of the block-diagonal matrix $A$.</span></div>
<div class="line"><a id="l00094" name="l00094"></a><span class="lineno"> 94</span><span class="comment"> ///</span></div>
<div class="line"><a id="l00095" name="l00095"></a><span class="lineno"> 95</span><span class="comment"> /// The LU factors and pivots of $A$ should be obtained by first calling</span></div>
<div class="line"><a id="l00096" name="l00096"></a><span class="lineno"> 96</span><span class="comment"> /// LUFactor(). $A$ has shape (n, n, n_mat) and $x$ has shape (n, n_rhs,</span></div>
<div class="line"><a id="l00097" name="l00097"></a><span class="lineno"> 97</span><span class="comment"> /// n_mat).</span></div>
<div class="line"><a id="l00098" name="l00098"></a><span class="lineno"> 98</span><span class="comment"> ///</span></div>
<div class="line"><a id="l00099" name="l00099"></a><span class="lineno"> 99</span><span class="comment"> /// @warning LUSolve() and LUFactor() should be called using the same backend</span></div>
<div class="line"><a id="l00100" name="l00100"></a><span class="lineno"> 100</span><span class="comment"> /// because of potential incompatibilities (e.g. 0-based or 1-based</span></div>
<div class="line"><a id="l00101" name="l00101"></a><span class="lineno"> 101</span><span class="comment"> /// indexing).</span></div>
<div class="line"><a id="l00102" name="l00102"></a><span class="lineno"> 102</span><span class="comment"></span> <span class="keyword">static</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classmfem_1_1BatchedLinAlg.html#adecb845877975afb36eba518c5cda331">LUSolve</a>(<span class="keyword">const</span> <a class="code hl_class" href="classmfem_1_1DenseTensor.html">DenseTensor</a> &A, <span class="keyword">const</span> <a class="code hl_class" href="classmfem_1_1Array.html">Array<int></a> &P, <a class="code hl_class" href="classmfem_1_1Vector.html">Vector</a> &x);<span class="comment"></span></div>
<div class="line"><a id="l00103" name="l00103"></a><span class="lineno"> 103</span><span class="comment"> /// @brief Returns true if the requested backend is available.</span></div>
<div class="line"><a id="l00104" name="l00104"></a><span class="lineno"> 104</span><span class="comment"> ///</span></div>
<div class="line"><a id="l00105" name="l00105"></a><span class="lineno"> 105</span><span class="comment"> /// The available backends depend on which third-party libraries MFEM is</span></div>
<div class="line"><a id="l00106" name="l00106"></a><span class="lineno"> 106</span><span class="comment"> /// compiled with, and whether the CUDA/HIP device is enabled.</span></div>
<div class="line"><a id="l00107" name="l00107"></a><span class="lineno"> 107</span><span class="comment"></span> <span class="keyword">static</span> <span class="keywordtype">bool</span> <a class="code hl_function" href="classmfem_1_1BatchedLinAlg.html#aa58c66de987a2bdcd8bfcfa37be585dd">IsAvailable</a>(<a class="code hl_struct" href="structmfem_1_1Backend.html">Backend</a> backend);<span class="comment"></span></div>
<div class="line"><a id="l00108" name="l00108"></a><span class="lineno"> 108</span><span class="comment"> /// Set the default backend for batched linear algebra operations.</span></div>
<div class="line"><a id="l00109" name="l00109"></a><span class="lineno"> 109</span><span class="comment"></span> <span class="keyword">static</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classmfem_1_1BatchedLinAlg.html#ad366fa4e973dae68161e151488c29303">SetActiveBackend</a>(<a class="code hl_struct" href="structmfem_1_1Backend.html">Backend</a> backend);<span class="comment"></span></div>
<div class="line"><a id="l00110" name="l00110"></a><span class="lineno"> 110</span><span class="comment"> /// Get the default backend for batched linear algebra operations.</span></div>
<div class="line"><a id="l00111" name="l00111"></a><span class="lineno"> 111</span><span class="comment"></span> <span class="keyword">static</span> <a class="code hl_struct" href="structmfem_1_1Backend.html">Backend</a> <a class="code hl_function" href="classmfem_1_1BatchedLinAlg.html#a9f37b31769ced7171dc6518ba814a009">GetActiveBackend</a>();<span class="comment"></span></div>
<div class="line"><a id="l00112" name="l00112"></a><span class="lineno"> 112</span><span class="comment"> /// @brief Get the BatchedLinAlgBase object associated with a specific</span></div>
<div class="line"><a id="l00113" name="l00113"></a><span class="lineno"> 113</span><span class="comment"> /// backend.</span></div>
<div class="line"><a id="l00114" name="l00114"></a><span class="lineno"> 114</span><span class="comment"> ///</span></div>
<div class="line"><a id="l00115" name="l00115"></a><span class="lineno"> 115</span><span class="comment"> /// This allows the user to perform specific operations with a backend</span></div>
<div class="line"><a id="l00116" name="l00116"></a><span class="lineno"> 116</span><span class="comment"> /// different from the active backend.</span></div>
<div class="line"><a id="l00117" name="l00117"></a><span class="lineno"> 117</span><span class="comment"></span> <span class="keyword">static</span> <span class="keyword">const</span> <a class="code hl_class" href="classmfem_1_1BatchedLinAlgBase.html">BatchedLinAlgBase</a> &<a class="code hl_function" href="classmfem_1_1BatchedLinAlg.html#a78a5b493f10cef0671e69add71fe4ae5">Get</a>(<a class="code hl_struct" href="structmfem_1_1Backend.html">Backend</a> backend);</div>
<div class="line"><a id="l00118" name="l00118"></a><span class="lineno"> 118</span>};</div>
</div>
<div class="line"><a id="l00119" name="l00119"></a><span class="lineno"> 119</span><span class="comment"></span> </div>
<div class="line"><a id="l00120" name="l00120"></a><span class="lineno"> 120</span><span class="comment">/// Abstract base clase for batched linear algebra operations.</span></div>
<div class="foldopen" id="foldopen00121" data-start="{" data-end="};">
<div class="line"><a id="l00121" name="l00121"></a><span class="lineno"><a class="line" href="classmfem_1_1BatchedLinAlgBase.html"> 121</a></span><span class="comment"></span><span class="keyword">class </span><a class="code hl_class" href="classmfem_1_1BatchedLinAlgBase.html">BatchedLinAlgBase</a></div>
<div class="line"><a id="l00122" name="l00122"></a><span class="lineno"> 122</span>{</div>
<div class="line"><a id="l00123" name="l00123"></a><span class="lineno"> 123</span><span class="keyword">public</span>:</div>
<div class="line"><a id="l00124" name="l00124"></a><span class="lineno"><a class="line" href="classmfem_1_1BatchedLinAlgBase.html#af8476c6d83feeb9c698b3b231ebb3506"> 124</a></span> <span class="keyword">using </span><a class="code hl_enumeration" href="classmfem_1_1BatchedLinAlg.html#ade3f74dd4281b3e06ab1cee7448b0d8d">Op</a> = <a class="code hl_enumeration" href="classmfem_1_1BatchedLinAlg.html#ade3f74dd4281b3e06ab1cee7448b0d8d">BatchedLinAlg::Op</a>;<span class="comment"></span></div>
<div class="line"><a id="l00125" name="l00125"></a><span class="lineno"> 125</span><span class="comment"> /// See BatchedLinAlg::AddMult.</span></div>
<div class="line"><a id="l00126" name="l00126"></a><span class="lineno"><a class="line" href="classmfem_1_1BatchedLinAlgBase.html#ae48ec53078353fe9ff41a40f5eef31bd"> 126</a></span><span class="comment"></span> <span class="keyword">virtual</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classmfem_1_1BatchedLinAlgBase.html#ae48ec53078353fe9ff41a40f5eef31bd">AddMult</a>(<span class="keyword">const</span> <a class="code hl_class" href="classmfem_1_1DenseTensor.html">DenseTensor</a> &A, <span class="keyword">const</span> <a class="code hl_class" href="classmfem_1_1Vector.html">Vector</a> &x, <a class="code hl_class" href="classmfem_1_1Vector.html">Vector</a> &y,</div>
<div class="line"><a id="l00127" name="l00127"></a><span class="lineno"> 127</span> <a class="code hl_typedef" href="namespacemfem.html#ac8ba85124b5de90bf368176b69d6019c">real_t</a> <a class="code hl_variable" href="ex15_8cpp.html#a7dc6cd505f60a008e9246be53a4d460d">alpha</a> = 1.0, <a class="code hl_typedef" href="namespacemfem.html#ac8ba85124b5de90bf368176b69d6019c">real_t</a> <a class="code hl_variable" href="convection-diffusion_8cpp.html#a20d35965a489ee0e4396d34b3867aaed">beta</a> = 1.0,</div>
<div class="line"><a id="l00128" name="l00128"></a><span class="lineno"> 128</span> <a class="code hl_enumeration" href="classmfem_1_1BatchedLinAlg.html#ade3f74dd4281b3e06ab1cee7448b0d8d">Op</a> op = Op::N) <span class="keyword">const</span> = 0;<span class="comment"></span></div>
<div class="line"><a id="l00129" name="l00129"></a><span class="lineno"> 129</span><span class="comment"> /// See BatchedLinAlg::Mult.</span></div>
<div class="line"><a id="l00130" name="l00130"></a><span class="lineno"> 130</span><span class="comment"></span> <span class="keyword">virtual</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classmfem_1_1BatchedLinAlgBase.html#a4c221f763b91b03dddca02221b7c69ac">Mult</a>(<span class="keyword">const</span> <a class="code hl_class" href="classmfem_1_1DenseTensor.html">DenseTensor</a> &A, <span class="keyword">const</span> <a class="code hl_class" href="classmfem_1_1Vector.html">Vector</a> &x, <a class="code hl_class" href="classmfem_1_1Vector.html">Vector</a> &y) <span class="keyword">const</span>;<span class="comment"></span></div>
<div class="line"><a id="l00131" name="l00131"></a><span class="lineno"> 131</span><span class="comment"> /// See BatchedLinAlg::MultTranspose.</span></div>
<div class="line"><a id="l00132" name="l00132"></a><span class="lineno"> 132</span><span class="comment"></span> <span class="keyword">virtual</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classmfem_1_1BatchedLinAlgBase.html#aa0fa69edca49f25f69ba57c223625568">MultTranspose</a>(<span class="keyword">const</span> <a class="code hl_class" href="classmfem_1_1DenseTensor.html">DenseTensor</a> &A, <span class="keyword">const</span> <a class="code hl_class" href="classmfem_1_1Vector.html">Vector</a> &x,</div>
<div class="line"><a id="l00133" name="l00133"></a><span class="lineno"> 133</span> <a class="code hl_class" href="classmfem_1_1Vector.html">Vector</a> &y) <span class="keyword">const</span>;<span class="comment"></span></div>
<div class="line"><a id="l00134" name="l00134"></a><span class="lineno"> 134</span><span class="comment"> /// See BatchedLinAlg::Invert.</span></div>
<div class="line"><a id="l00135" name="l00135"></a><span class="lineno"><a class="line" href="classmfem_1_1BatchedLinAlgBase.html#a68bff05aa750628e1a38b8a6e5cd9fbe"> 135</a></span><span class="comment"></span> <span class="keyword">virtual</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classmfem_1_1BatchedLinAlgBase.html#a68bff05aa750628e1a38b8a6e5cd9fbe">Invert</a>(<a class="code hl_class" href="classmfem_1_1DenseTensor.html">DenseTensor</a> &A) <span class="keyword">const</span> = 0;<span class="comment"></span></div>
<div class="line"><a id="l00136" name="l00136"></a><span class="lineno"> 136</span><span class="comment"> /// See BatchedLinAlg::LUFactor.</span></div>
<div class="line"><a id="l00137" name="l00137"></a><span class="lineno"><a class="line" href="classmfem_1_1BatchedLinAlgBase.html#a6bec58c978e5be6eea04d5b9ad0fd917"> 137</a></span><span class="comment"></span> <span class="keyword">virtual</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classmfem_1_1BatchedLinAlgBase.html#a6bec58c978e5be6eea04d5b9ad0fd917">LUFactor</a>(<a class="code hl_class" href="classmfem_1_1DenseTensor.html">DenseTensor</a> &A, <a class="code hl_class" href="classmfem_1_1Array.html">Array<int></a> &P) <span class="keyword">const</span> = 0;<span class="comment"></span></div>
<div class="line"><a id="l00138" name="l00138"></a><span class="lineno"> 138</span><span class="comment"> /// See BatchedLinAlg::LUSolve.</span></div>
<div class="line"><a id="l00139" name="l00139"></a><span class="lineno"><a class="line" href="classmfem_1_1BatchedLinAlgBase.html#acc73422d13a2b0da92f8830514a3fa7e"> 139</a></span><span class="comment"></span> <span class="keyword">virtual</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classmfem_1_1BatchedLinAlgBase.html#acc73422d13a2b0da92f8830514a3fa7e">LUSolve</a>(<span class="keyword">const</span> <a class="code hl_class" href="classmfem_1_1DenseTensor.html">DenseTensor</a> &LU, <span class="keyword">const</span> <a class="code hl_class" href="classmfem_1_1Array.html">Array<int></a> &P,</div>
<div class="line"><a id="l00140" name="l00140"></a><span class="lineno"> 140</span> <a class="code hl_class" href="classmfem_1_1Vector.html">Vector</a> &x) <span class="keyword">const</span> = 0;<span class="comment"></span></div>
<div class="line"><a id="l00141" name="l00141"></a><span class="lineno"> 141</span><span class="comment"> /// Virtual destructor.</span></div>
<div class="line"><a id="l00142" name="l00142"></a><span class="lineno"><a class="line" href="classmfem_1_1BatchedLinAlgBase.html#a5bb3598844830353c74b887e607b2419"> 142</a></span><span class="comment"></span> <span class="keyword">virtual</span> <a class="code hl_function" href="classmfem_1_1BatchedLinAlgBase.html#a5bb3598844830353c74b887e607b2419">~BatchedLinAlgBase</a>() { }</div>
<div class="line"><a id="l00143" name="l00143"></a><span class="lineno"> 143</span>};</div>
</div>
<div class="line"><a id="l00144" name="l00144"></a><span class="lineno"> 144</span> </div>
<div class="line"><a id="l00145" name="l00145"></a><span class="lineno"> 145</span>} <span class="comment">// namespace mfem</span></div>
<div class="line"><a id="l00146" name="l00146"></a><span class="lineno"> 146</span> </div>
<div class="line"><a id="l00147" name="l00147"></a><span class="lineno"> 147</span><span class="preprocessor">#endif</span></div>
<div class="ttc" id="aclassmfem_1_1Array_html"><div class="ttname"><a href="classmfem_1_1Array.html">mfem::Array</a></div><div class="ttdef"><b>Definition</b> <a href="array_8hpp_source.html#l00046">array.hpp:47</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlgBase_html"><div class="ttname"><a href="classmfem_1_1BatchedLinAlgBase.html">mfem::BatchedLinAlgBase</a></div><div class="ttdoc">Abstract base clase for batched linear algebra operations.</div><div class="ttdef"><b>Definition</b> <a href="#l00121">batched.hpp:122</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlgBase_html_a4c221f763b91b03dddca02221b7c69ac"><div class="ttname"><a href="classmfem_1_1BatchedLinAlgBase.html#a4c221f763b91b03dddca02221b7c69ac">mfem::BatchedLinAlgBase::Mult</a></div><div class="ttdeci">virtual void Mult(const DenseTensor &A, const Vector &x, Vector &y) const</div><div class="ttdoc">See BatchedLinAlg::Mult.</div><div class="ttdef"><b>Definition</b> <a href="batched_8cpp_source.html#l00110">batched.cpp:110</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlgBase_html_a5bb3598844830353c74b887e607b2419"><div class="ttname"><a href="classmfem_1_1BatchedLinAlgBase.html#a5bb3598844830353c74b887e607b2419">mfem::BatchedLinAlgBase::~BatchedLinAlgBase</a></div><div class="ttdeci">virtual ~BatchedLinAlgBase()</div><div class="ttdoc">Virtual destructor.</div><div class="ttdef"><b>Definition</b> <a href="#l00142">batched.hpp:142</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlgBase_html_a68bff05aa750628e1a38b8a6e5cd9fbe"><div class="ttname"><a href="classmfem_1_1BatchedLinAlgBase.html#a68bff05aa750628e1a38b8a6e5cd9fbe">mfem::BatchedLinAlgBase::Invert</a></div><div class="ttdeci">virtual void Invert(DenseTensor &A) const =0</div><div class="ttdoc">See BatchedLinAlg::Invert.</div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlgBase_html_a6bec58c978e5be6eea04d5b9ad0fd917"><div class="ttname"><a href="classmfem_1_1BatchedLinAlgBase.html#a6bec58c978e5be6eea04d5b9ad0fd917">mfem::BatchedLinAlgBase::LUFactor</a></div><div class="ttdeci">virtual void LUFactor(DenseTensor &A, Array< int > &P) const =0</div><div class="ttdoc">See BatchedLinAlg::LUFactor.</div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlgBase_html_aa0fa69edca49f25f69ba57c223625568"><div class="ttname"><a href="classmfem_1_1BatchedLinAlgBase.html#aa0fa69edca49f25f69ba57c223625568">mfem::BatchedLinAlgBase::MultTranspose</a></div><div class="ttdeci">virtual void MultTranspose(const DenseTensor &A, const Vector &x, Vector &y) const</div><div class="ttdoc">See BatchedLinAlg::MultTranspose.</div><div class="ttdef"><b>Definition</b> <a href="batched_8cpp_source.html#l00116">batched.cpp:116</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlgBase_html_acc73422d13a2b0da92f8830514a3fa7e"><div class="ttname"><a href="classmfem_1_1BatchedLinAlgBase.html#acc73422d13a2b0da92f8830514a3fa7e">mfem::BatchedLinAlgBase::LUSolve</a></div><div class="ttdeci">virtual void LUSolve(const DenseTensor &LU, const Array< int > &P, Vector &x) const =0</div><div class="ttdoc">See BatchedLinAlg::LUSolve.</div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlgBase_html_ae48ec53078353fe9ff41a40f5eef31bd"><div class="ttname"><a href="classmfem_1_1BatchedLinAlgBase.html#ae48ec53078353fe9ff41a40f5eef31bd">mfem::BatchedLinAlgBase::AddMult</a></div><div class="ttdeci">virtual void AddMult(const DenseTensor &A, const Vector &x, Vector &y, real_t alpha=1.0, real_t beta=1.0, Op op=Op::N) const =0</div><div class="ttdoc">See BatchedLinAlg::AddMult.</div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html">mfem::BatchedLinAlg</a></div><div class="ttdoc">Class for performing batched linear algebra operations, potentially using accelerated algorithms (GPU...</div><div class="ttdef"><b>Definition</b> <a href="#l00031">batched.hpp:32</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_a327ee59ae1496af66d18207a94d491fd"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#a327ee59ae1496af66d18207a94d491fd">mfem::BatchedLinAlg::Mult</a></div><div class="ttdeci">static void Mult(const DenseTensor &A, const Vector &x, Vector &y)</div><div class="ttdoc">Computes (e.g. by calling AddMult(A,x,y,1,0,Op::N)).</div><div class="ttdef"><b>Definition</b> <a href="batched_8cpp_source.html#l00060">batched.cpp:60</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_a40da089d2e89a7c7ecef22423473b6ca"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#a40da089d2e89a7c7ecef22423473b6ca">mfem::BatchedLinAlg::Backend</a></div><div class="ttdeci">Backend</div><div class="ttdoc">Available backends for implementations of batched algorithms.</div><div class="ttdef"><b>Definition</b> <a href="#l00038">batched.hpp:39</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_a40da089d2e89a7c7ecef22423473b6caa388e52e74ac0158131847fcc330b16ee"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#a40da089d2e89a7c7ecef22423473b6caa388e52e74ac0158131847fcc330b16ee">mfem::BatchedLinAlg::GPU_BLAS</a></div><div class="ttdeci">@ GPU_BLAS</div><div class="ttdoc">Either cuBLAS or hipBLAS, depending on whether MFEM is using CUDA or HIP. Not available otherwise.</div><div class="ttdef"><b>Definition</b> <a href="#l00045">batched.hpp:45</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_a40da089d2e89a7c7ecef22423473b6caa462cb27cf8966b995bce9663e7d327ce"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#a40da089d2e89a7c7ecef22423473b6caa462cb27cf8966b995bce9663e7d327ce">mfem::BatchedLinAlg::MAGMA</a></div><div class="ttdeci">@ MAGMA</div><div class="ttdoc">MAGMA backend, only available if MFEM is compiled with MAGMA support.</div><div class="ttdef"><b>Definition</b> <a href="#l00047">batched.hpp:47</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_a40da089d2e89a7c7ecef22423473b6caa467cc6ba51d7366f915ac4413ab5ef09"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#a40da089d2e89a7c7ecef22423473b6caa467cc6ba51d7366f915ac4413ab5ef09">mfem::BatchedLinAlg::NUM_BACKENDS</a></div><div class="ttdeci">@ NUM_BACKENDS</div><div class="ttdoc">Counter for the number of backends.</div><div class="ttdef"><b>Definition</b> <a href="#l00050">batched.hpp:49</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_a40da089d2e89a7c7ecef22423473b6caadd4da6e88f75868880f6692ed6b5eeb4"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#a40da089d2e89a7c7ecef22423473b6caadd4da6e88f75868880f6692ed6b5eeb4">mfem::BatchedLinAlg::NATIVE</a></div><div class="ttdeci">@ NATIVE</div><div class="ttdoc">The standard MFEM backend, implemented using mfem::forall kernels. Not as performant as the other ker...</div><div class="ttdef"><b>Definition</b> <a href="#l00042">batched.hpp:42</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_a78a5b493f10cef0671e69add71fe4ae5"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#a78a5b493f10cef0671e69add71fe4ae5">mfem::BatchedLinAlg::Get</a></div><div class="ttdeci">static const BatchedLinAlgBase & Get(Backend backend)</div><div class="ttdoc">Get the BatchedLinAlgBase object associated with a specific backend.</div><div class="ttdef"><b>Definition</b> <a href="batched_8cpp_source.html#l00103">batched.cpp:103</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_a806592aaaa535f422453f75e855735ea"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#a806592aaaa535f422453f75e855735ea">mfem::BatchedLinAlg::MultTranspose</a></div><div class="ttdeci">static void MultTranspose(const DenseTensor &A, const Vector &x, Vector &y)</div><div class="ttdoc">Computes (e.g. by calling AddMult(A,x,y,1,0,Op::T)).</div><div class="ttdef"><b>Definition</b> <a href="batched_8cpp_source.html#l00065">batched.cpp:65</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_a9f37b31769ced7171dc6518ba814a009"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#a9f37b31769ced7171dc6518ba814a009">mfem::BatchedLinAlg::GetActiveBackend</a></div><div class="ttdeci">static Backend GetActiveBackend()</div><div class="ttdoc">Get the default backend for batched linear algebra operations.</div><div class="ttdef"><b>Definition</b> <a href="batched_8cpp_source.html#l00098">batched.cpp:98</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_aa58c66de987a2bdcd8bfcfa37be585dd"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#aa58c66de987a2bdcd8bfcfa37be585dd">mfem::BatchedLinAlg::IsAvailable</a></div><div class="ttdeci">static bool IsAvailable(Backend backend)</div><div class="ttdoc">Returns true if the requested backend is available.</div><div class="ttdef"><b>Definition</b> <a href="batched_8cpp_source.html#l00087">batched.cpp:87</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_aaf20d258a08692f4e9c7e3125d980658"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#aaf20d258a08692f4e9c7e3125d980658">mfem::BatchedLinAlg::LUFactor</a></div><div class="ttdeci">static void LUFactor(DenseTensor &A, Array< int > &P)</div><div class="ttdoc">Replaces the block diagonal matrix with its LU factors. The pivots are stored in P.</div><div class="ttdef"><b>Definition</b> <a href="batched_8cpp_source.html#l00076">batched.cpp:76</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_ab3676ca00da6882e00dfb7718b5241dd"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#ab3676ca00da6882e00dfb7718b5241dd">mfem::BatchedLinAlg::AddMult</a></div><div class="ttdeci">static void AddMult(const DenseTensor &A, const Vector &x, Vector &y, real_t alpha=1.0, real_t beta=1.0, Op op=Op::N)</div><div class="ttdoc">Computes .</div><div class="ttdef"><b>Definition</b> <a href="batched_8cpp_source.html#l00054">batched.cpp:54</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_ad366fa4e973dae68161e151488c29303"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#ad366fa4e973dae68161e151488c29303">mfem::BatchedLinAlg::SetActiveBackend</a></div><div class="ttdeci">static void SetActiveBackend(Backend backend)</div><div class="ttdoc">Set the default backend for batched linear algebra operations.</div><div class="ttdef"><b>Definition</b> <a href="batched_8cpp_source.html#l00092">batched.cpp:92</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_ade3f74dd4281b3e06ab1cee7448b0d8d"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#ade3f74dd4281b3e06ab1cee7448b0d8d">mfem::BatchedLinAlg::Op</a></div><div class="ttdeci">Op</div><div class="ttdoc">Operation type (transposed or not transposed)</div><div class="ttdef"><b>Definition</b> <a href="#l00053">batched.hpp:54</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_ade3f74dd4281b3e06ab1cee7448b0d8da19710f7a1d92fa9659d19cef79fe781b"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#ade3f74dd4281b3e06ab1cee7448b0d8da19710f7a1d92fa9659d19cef79fe781b">mfem::BatchedLinAlg::T</a></div><div class="ttdeci">@ T</div><div class="ttdoc">Transposed.</div><div class="ttdef"><b>Definition</b> <a href="#l00057">batched.hpp:56</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_ade3f74dd4281b3e06ab1cee7448b0d8dabefd6c3489b5cd88b472b82cff88df8f"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#ade3f74dd4281b3e06ab1cee7448b0d8dabefd6c3489b5cd88b472b82cff88df8f">mfem::BatchedLinAlg::N</a></div><div class="ttdeci">@ N</div><div class="ttdoc">Not transposed.</div><div class="ttdef"><b>Definition</b> <a href="#l00055">batched.hpp:55</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_adecb845877975afb36eba518c5cda331"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#adecb845877975afb36eba518c5cda331">mfem::BatchedLinAlg::LUSolve</a></div><div class="ttdeci">static void LUSolve(const DenseTensor &A, const Array< int > &P, Vector &x)</div><div class="ttdoc">Replaces with , given the LU factors A and pivots P of the block-diagonal matrix .</div><div class="ttdef"><b>Definition</b> <a href="batched_8cpp_source.html#l00081">batched.cpp:81</a></div></div>
<div class="ttc" id="aclassmfem_1_1BatchedLinAlg_html_ae64e57da2451c570f6e0fa06ebd8a111"><div class="ttname"><a href="classmfem_1_1BatchedLinAlg.html#ae64e57da2451c570f6e0fa06ebd8a111">mfem::BatchedLinAlg::Invert</a></div><div class="ttdeci">static void Invert(DenseTensor &A)</div><div class="ttdoc">Replaces the block diagonal matrix with its inverse .</div><div class="ttdef"><b>Definition</b> <a href="batched_8cpp_source.html#l00071">batched.cpp:71</a></div></div>
<div class="ttc" id="aclassmfem_1_1DenseTensor_html"><div class="ttname"><a href="classmfem_1_1DenseTensor.html">mfem::DenseTensor</a></div><div class="ttdoc">Rank 3 tensor (array of matrices)</div><div class="ttdef"><b>Definition</b> <a href="densemat_8hpp_source.html#l01111">densemat.hpp:1112</a></div></div>
<div class="ttc" id="aclassmfem_1_1Vector_html"><div class="ttname"><a href="classmfem_1_1Vector.html">mfem::Vector</a></div><div class="ttdoc">Vector data type.</div><div class="ttdef"><b>Definition</b> <a href="vector_8hpp_source.html#l00081">vector.hpp:82</a></div></div>
<div class="ttc" id="aconfig_8hpp_html"><div class="ttname"><a href="config_8hpp.html">config.hpp</a></div></div>
<div class="ttc" id="aconvection-diffusion_8cpp_html_a20d35965a489ee0e4396d34b3867aaed"><div class="ttname"><a href="convection-diffusion_8cpp.html#a20d35965a489ee0e4396d34b3867aaed">beta</a></div><div class="ttdeci">Vector beta</div><div class="ttdef"><b>Definition</b> <a href="convection-diffusion_8cpp_source.html#l00082">convection-diffusion.cpp:82</a></div></div>
<div class="ttc" id="adensemat_8hpp_html"><div class="ttname"><a href="densemat_8hpp.html">densemat.hpp</a></div></div>
<div class="ttc" id="aex15_8cpp_html_a7dc6cd505f60a008e9246be53a4d460d"><div class="ttname"><a href="ex15_8cpp.html#a7dc6cd505f60a008e9246be53a4d460d">alpha</a></div><div class="ttdeci">const real_t alpha</div><div class="ttdef"><b>Definition</b> <a href="ex15_8cpp_source.html#l00369">ex15.cpp:369</a></div></div>
<div class="ttc" id="anamespacemfem_html"><div class="ttname"><a href="namespacemfem.html">mfem</a></div><div class="ttdef"><b>Definition</b> <a href="CodeDocumentation_8dox_source.html#l00001">CodeDocumentation.dox:1</a></div></div>
<div class="ttc" id="anamespacemfem_html_ac8ba85124b5de90bf368176b69d6019c"><div class="ttname"><a href="namespacemfem.html#ac8ba85124b5de90bf368176b69d6019c">mfem::real_t</a></div><div class="ttdeci">float real_t</div><div class="ttdef"><b>Definition</b> <a href="config_8hpp_source.html#l00043">config.hpp:43</a></div></div>
<div class="ttc" id="astructmfem_1_1Backend_html"><div class="ttname"><a href="structmfem_1_1Backend.html">mfem::Backend</a></div><div class="ttdoc">MFEM backends.</div><div class="ttdef"><b>Definition</b> <a href="device_8hpp_source.html#l00027">device.hpp:28</a></div></div>
</div><!-- fragment --></div><!-- contents -->
<!-- start footer part -->
<hr class="footer"/><address class="footer"><small>
Generated on Wed Apr 9 2025 16:53:04 for MFEM by <a href="https://www.doxygen.org/index.html"><img class="footer" src="doxygen.svg" width="104" height="31" alt="doxygen"/></a> 1.11.0
</small></address>
</div><!-- doc-content -->
</body>
</html>