OpenFOAM: API Guide
v2006
The open source CFD toolbox
fvcLaplacian.C
Go to the documentation of this file.
1
/*---------------------------------------------------------------------------*\
2
========= |
3
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4
\\ / O peration |
5
\\ / A nd | www.openfoam.com
6
\\/ M anipulation |
7
-------------------------------------------------------------------------------
8
Copyright (C) 2011-2016 OpenFOAM Foundation
9
-------------------------------------------------------------------------------
10
License
11
This file is part of OpenFOAM.
12
13
OpenFOAM is free software: you can redistribute it and/or modify it
14
under the terms of the GNU General Public License as published by
15
the Free Software Foundation, either version 3 of the License, or
16
(at your option) any later version.
17
18
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21
for more details.
22
23
You should have received a copy of the GNU General Public License
24
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26
\*---------------------------------------------------------------------------*/
27
28
#include "
fvcLaplacian.H
"
29
#include "
fvMesh.H
"
30
#include "
laplacianScheme.H
"
31
32
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34
namespace
Foam
35
{
36
37
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
39
namespace
fvc
40
{
41
42
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43
44
template
<
class
Type>
45
tmp<GeometricField<Type, fvPatchField, volMesh>>
46
laplacian
47
(
48
const
GeometricField<Type, fvPatchField, volMesh>
& vf,
49
const
word
&
name
50
)
51
{
52
return
fv::laplacianScheme<Type, scalar>::New
53
(
54
vf.mesh(),
55
vf.mesh().laplacianScheme(
name
)
56
).
ref
().fvcLaplacian(vf);
57
}
58
59
60
template
<
class
Type>
61
tmp<GeometricField<Type, fvPatchField, volMesh>
>
62
laplacian
63
(
64
const
tmp
<
GeometricField<Type, fvPatchField, volMesh>
>& tvf,
65
const
word
&
name
66
)
67
{
68
tmp<GeometricField<Type, fvPatchField, volMesh>
> Laplacian
69
(
70
fvc::laplacian
(tvf(),
name
)
71
);
72
tvf.clear();
73
return
Laplacian;
74
}
75
76
77
template
<
class
Type>
78
tmp<GeometricField<Type, fvPatchField, volMesh>
>
79
laplacian
80
(
81
const
GeometricField<Type, fvPatchField, volMesh>
& vf
82
)
83
{
84
return
fvc::laplacian
(vf,
"laplacian("
+ vf.name() +
')'
);
85
}
86
87
88
template
<
class
Type>
89
tmp<GeometricField<Type, fvPatchField, volMesh>
>
90
laplacian
91
(
92
const
tmp
<
GeometricField<Type, fvPatchField, volMesh>
>& tvf
93
)
94
{
95
tmp<GeometricField<Type, fvPatchField, volMesh>
> Laplacian
96
(
97
fvc::laplacian
(tvf())
98
);
99
tvf.clear();
100
return
Laplacian;
101
}
102
103
104
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
105
106
template
<
class
Type,
class
GType>
107
tmp<GeometricField<Type, fvPatchField, volMesh>
>
108
laplacian
109
(
110
const
dimensioned<GType>
&
gamma
,
111
const
GeometricField<Type, fvPatchField, volMesh>
& vf,
112
const
word
&
name
113
)
114
{
115
GeometricField<GType, fvsPatchField, surfaceMesh>
Gamma
116
(
117
IOobject
118
(
119
gamma
.name(),
120
vf.instance(),
121
vf.mesh(),
122
IOobject::NO_READ
123
),
124
vf.mesh(),
125
gamma
126
);
127
128
return
fvc::laplacian
(Gamma, vf,
name
);
129
}
130
131
132
template
<
class
Type,
class
GType>
133
tmp<GeometricField<Type, fvPatchField, volMesh>
>
134
laplacian
135
(
136
const
dimensioned<GType>
&
gamma
,
137
const
tmp
<
GeometricField<Type, fvPatchField, volMesh>
>& tvf,
138
const
word
&
name
139
)
140
{
141
tmp<GeometricField<Type, fvPatchField, volMesh>
> Laplacian
142
(
143
fvc::laplacian
(
gamma
, tvf(),
name
)
144
);
145
tvf.clear();
146
return
Laplacian;
147
}
148
149
150
template
<
class
Type,
class
GType>
151
tmp<GeometricField<Type, fvPatchField, volMesh>
>
152
laplacian
153
(
154
const
dimensioned<GType>
&
gamma
,
155
const
GeometricField<Type, fvPatchField, volMesh>
& vf
156
)
157
{
158
GeometricField<GType, fvsPatchField, surfaceMesh>
Gamma
159
(
160
IOobject
161
(
162
gamma
.name(),
163
vf.instance(),
164
vf.mesh(),
165
IOobject::NO_READ
166
),
167
vf.mesh(),
168
gamma
169
);
170
171
return
fvc::laplacian
(Gamma, vf);
172
}
173
174
175
template
<
class
Type,
class
GType>
176
tmp<GeometricField<Type, fvPatchField, volMesh>
>
177
laplacian
178
(
179
const
dimensioned<GType>
&
gamma
,
180
const
tmp
<
GeometricField<Type, fvPatchField, volMesh>
>& tvf
181
)
182
{
183
tmp<GeometricField<Type, fvPatchField, volMesh>
> Laplacian
184
(
185
fvc::laplacian
(
gamma
, tvf())
186
);
187
tvf.clear();
188
return
Laplacian;
189
}
190
191
192
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193
194
template
<
class
Type,
class
GType>
195
tmp<GeometricField<Type, fvPatchField, volMesh>
>
196
laplacian
197
(
198
const
GeometricField<GType, fvPatchField, volMesh>
&
gamma
,
199
const
GeometricField<Type, fvPatchField, volMesh>
& vf,
200
const
word
&
name
201
)
202
{
203
return
fv::laplacianScheme<Type, GType>::New
204
(
205
vf.mesh(),
206
vf.mesh().laplacianScheme(
name
)
207
).
ref
().fvcLaplacian(
gamma
, vf);
208
}
209
210
211
template
<
class
Type,
class
GType>
212
tmp<GeometricField<Type, fvPatchField, volMesh>
>
213
laplacian
214
(
215
const
tmp
<
GeometricField<GType, fvPatchField, volMesh>
>& tgamma,
216
const
GeometricField<Type, fvPatchField, volMesh>
& vf,
217
const
word
&
name
218
)
219
{
220
tmp<GeometricField<Type, fvPatchField, volMesh>
> Laplacian
221
(
222
fvc::laplacian
(tgamma(), vf,
name
)
223
);
224
tgamma.clear();
225
return
Laplacian;
226
}
227
228
229
template
<
class
Type,
class
GType>
230
tmp<GeometricField<Type, fvPatchField, volMesh>
>
231
laplacian
232
(
233
const
GeometricField<GType, fvPatchField, volMesh>
&
gamma
,
234
const
tmp
<
GeometricField<Type, fvPatchField, volMesh>
>& tvf,
235
const
word
&
name
236
)
237
{
238
tmp<GeometricField<Type, fvPatchField, volMesh>
> Laplacian
239
(
240
fvc::laplacian
(
gamma
, tvf(),
name
)
241
);
242
tvf.clear();
243
return
Laplacian;
244
}
245
246
247
template
<
class
Type,
class
GType>
248
tmp<GeometricField<Type, fvPatchField, volMesh>
>
249
laplacian
250
(
251
const
tmp
<
GeometricField<GType, fvPatchField, volMesh>
>& tgamma,
252
const
tmp
<
GeometricField<Type, fvPatchField, volMesh>
>& tvf,
253
const
word
&
name
254
)
255
{
256
tmp<GeometricField<Type, fvPatchField, volMesh>
> Laplacian
257
(
258
fvc::laplacian
(tgamma(), tvf(),
name
)
259
);
260
tgamma.clear();
261
tvf.clear();
262
return
Laplacian;
263
}
264
265
266
template
<
class
Type,
class
GType>
267
tmp<GeometricField<Type, fvPatchField, volMesh>
>
268
laplacian
269
(
270
const
GeometricField<GType, fvPatchField, volMesh>
&
gamma
,
271
const
GeometricField<Type, fvPatchField, volMesh>
& vf
272
)
273
{
274
return
fvc::laplacian
275
(
276
gamma
,
277
vf,
278
"laplacian("
+
gamma
.name() +
','
+ vf.name() +
')'
279
);
280
}
281
282
283
template
<
class
Type,
class
GType>
284
tmp<GeometricField<Type, fvPatchField, volMesh>
>
285
laplacian
286
(
287
const
tmp
<
GeometricField<GType, fvPatchField, volMesh>
>& tgamma,
288
const
GeometricField<Type, fvPatchField, volMesh>
& vf
289
)
290
{
291
return
fvc::laplacian
292
(
293
tgamma,
294
vf,
295
"laplacian("
+ tgamma().
name
() +
','
+ vf.name() +
')'
296
);
297
}
298
299
300
template
<
class
Type,
class
GType>
301
tmp<GeometricField<Type, fvPatchField, volMesh>
>
302
laplacian
303
(
304
const
GeometricField<GType, fvPatchField, volMesh>
&
gamma
,
305
const
tmp
<
GeometricField<Type, fvPatchField, volMesh>
>& tvf
306
)
307
{
308
return
fvc::laplacian
309
(
310
gamma
,
311
tvf,
312
"laplacian("
+
gamma
.name() +
','
+ tvf().name() +
')'
313
);
314
}
315
316
317
template
<
class
Type,
class
GType>
318
tmp<GeometricField<Type, fvPatchField, volMesh>
>
319
laplacian
320
(
321
const
tmp
<
GeometricField<GType, fvPatchField, volMesh>
>& tgamma,
322
const
tmp
<
GeometricField<Type, fvPatchField, volMesh>
>& tvf
323
)
324
{
325
return
fvc::laplacian
326
(
327
tgamma,
328
tvf,
329
"laplacian("
+ tgamma().
name
() +
','
+ tvf().
name
() +
')'
330
);
331
}
332
333
334
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
335
336
template
<
class
Type,
class
GType>
337
tmp<GeometricField<Type, fvPatchField, volMesh>
>
338
laplacian
339
(
340
const
GeometricField<GType, fvsPatchField, surfaceMesh>
&
gamma
,
341
const
GeometricField<Type, fvPatchField, volMesh>
& vf,
342
const
word
&
name
343
)
344
{
345
return
fv::laplacianScheme<Type, GType>::New
346
(
347
vf.mesh(),
348
vf.mesh().laplacianScheme(
name
)
349
).
ref
().fvcLaplacian(
gamma
, vf);
350
}
351
352
353
template
<
class
Type,
class
GType>
354
tmp<GeometricField<Type, fvPatchField, volMesh>
>
355
laplacian
356
(
357
const
tmp
<
GeometricField<GType, fvsPatchField, surfaceMesh>
>& tgamma,
358
const
GeometricField<Type, fvPatchField, volMesh>
& vf,
359
const
word
&
name
360
)
361
{
362
tmp<GeometricField<Type, fvPatchField, volMesh>
> Laplacian
363
(
364
fvc::laplacian
(tgamma(), vf,
name
)
365
);
366
tgamma.clear();
367
return
Laplacian;
368
}
369
370
371
template
<
class
Type,
class
GType>
372
tmp<GeometricField<Type, fvPatchField, volMesh>
>
373
laplacian
374
(
375
const
GeometricField<GType, fvsPatchField, surfaceMesh>
&
gamma
,
376
const
tmp
<
GeometricField<Type, fvPatchField, volMesh>
>& tvf,
377
const
word
&
name
378
)
379
{
380
tmp<GeometricField<Type, fvPatchField, volMesh>
> Laplacian
381
(
382
fvc::laplacian
(
gamma
, tvf(),
name
)
383
);
384
tvf.clear();
385
return
Laplacian;
386
}
387
388
389
template
<
class
Type,
class
GType>
390
tmp<GeometricField<Type, fvPatchField, volMesh>
>
laplacian
391
(
392
const
tmp
<
GeometricField<GType, fvsPatchField, surfaceMesh>
>& tgamma,
393
const
tmp
<
GeometricField<Type, fvPatchField, volMesh>
>& tvf,
394
const
word
&
name
395
)
396
{
397
tmp<GeometricField<Type, fvPatchField, volMesh>
> Laplacian
398
(
399
fvc::laplacian
(tgamma(), tvf(),
name
)
400
);
401
tgamma.clear();
402
tvf.clear();
403
return
Laplacian;
404
}
405
406
407
template
<
class
Type,
class
GType>
408
tmp<GeometricField<Type, fvPatchField, volMesh>
>
409
laplacian
410
(
411
const
GeometricField<GType, fvsPatchField, surfaceMesh>
&
gamma
,
412
const
GeometricField<Type, fvPatchField, volMesh>
& vf
413
)
414
{
415
return
fvc::laplacian
416
(
417
gamma
,
418
vf,
419
"laplacian("
+
gamma
.name() +
','
+ vf.name() +
')'
420
);
421
}
422
423
424
template
<
class
Type,
class
GType>
425
tmp<GeometricField<Type, fvPatchField, volMesh>
>
426
laplacian
427
(
428
const
tmp
<
GeometricField<GType, fvsPatchField, surfaceMesh>
>& tgamma,
429
const
GeometricField<Type, fvPatchField, volMesh>
& vf
430
)
431
{
432
tmp<GeometricField<Type, fvPatchField, volMesh>
> Laplacian
433
(
434
fvc::laplacian
(tgamma(), vf)
435
);
436
tgamma.clear();
437
return
Laplacian;
438
}
439
440
441
template
<
class
Type,
class
GType>
442
tmp<GeometricField<Type, fvPatchField, volMesh>
>
443
laplacian
444
(
445
const
GeometricField<GType, fvsPatchField, surfaceMesh>
&
gamma
,
446
const
tmp
<
GeometricField<Type, fvPatchField, volMesh>
>& tvf
447
)
448
{
449
tmp<GeometricField<Type, fvPatchField, volMesh>
> Laplacian
450
(
451
fvc::laplacian
(
gamma
, tvf())
452
);
453
tvf.clear();
454
return
Laplacian;
455
}
456
457
458
template
<
class
Type,
class
GType>
459
tmp<GeometricField<Type, fvPatchField, volMesh>
>
laplacian
460
(
461
const
tmp
<
GeometricField<GType, fvsPatchField, surfaceMesh>
>& tgamma,
462
const
tmp
<
GeometricField<Type, fvPatchField, volMesh>
>& tvf
463
)
464
{
465
tmp<GeometricField<Type, fvPatchField, volMesh>
> Laplacian
466
(
467
fvc::laplacian
(tgamma(), tvf())
468
);
469
tgamma.clear();
470
tvf.clear();
471
return
Laplacian;
472
}
473
474
475
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
476
477
}
// End namespace fvc
478
479
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
480
481
}
// End namespace Foam
482
483
// ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition:
IOobject.H:104
Foam::fv::laplacianScheme::New
static tmp< laplacianScheme< Type, GType > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new laplacianScheme created on freestore.
Definition:
laplacianScheme.C:48
Foam::word
A class for handling words, derived from Foam::string.
Definition:
word.H:62
Foam::tmp
A class for managing temporary objects.
Definition:
PtrList.H:59
laplacianScheme.H
ref
rDeltaT ref()
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition:
fvcLaplacian.C:47
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition:
complex.C:76
Foam::dimensioned
Generic dimensioned Type class.
Definition:
dimensionedScalarFwd.H:43
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition:
atmBoundaryLayer.C:33
fvcLaplacian.H
Calculate the laplacian of the given field.
gamma
const scalar gamma
Definition:
EEqn.H:9
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition:
IOobject.H:123
src
finiteVolume
finiteVolume
fvc
fvcLaplacian.C
Generated by
1.8.17
OPENFOAM® is a registered
trademark
of OpenCFD Ltd.