facLaplacian.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) 2016-2017 Wikki Ltd
9-------------------------------------------------------------------------------
10License
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 "facLaplacian.H"
29#include "faMesh.H"
30#include "faLaplacianScheme.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
36
37// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
39namespace fac
40{
41
42// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43
44template<class Type>
45tmp<GeometricField<Type, faPatchField, areaMesh>>
47(
49 const word& name
50)
51{
53 (
54 vf.mesh(),
55 vf.mesh().laplacianScheme(name)
56 ).ref().facLaplacian(vf);
57}
58
59template<class Type>
62(
64 const word& name
65)
66{
68 (
69 fac::laplacian(tvf(), name)
70 );
71 tvf.clear();
72 return Laplacian;
73}
74
75
76template<class Type>
79(
81)
82{
83 return fac::laplacian(vf, "laplacian(" + vf.name() + ')');
84}
85
86template<class Type>
89(
91)
92{
94 (
95 fac::laplacian(tvf())
96 );
97 tvf.clear();
98 return Laplacian;
99}
100
101
102template<class Type>
105(
108 const word& name
109)
110{
111 return gamma*fac::laplacian(vf, name);
112}
113
114template<class Type>
117(
120 const word& name
121)
122{
124 (
125 fac::laplacian(gamma, tvf(), name)
126 );
127 tvf.clear();
128 return Laplacian;
129}
130
131
132template<class Type>
135(
138)
139{
140 return gamma*fac::laplacian
141 (
142 vf, "laplacian(" + gamma.name() + ',' + vf.name() + ')'
143 );
144}
145
146template<class Type>
149(
152)
153{
155 (
156 fac::laplacian(gamma, tvf())
157 );
158 tvf.clear();
159 return Laplacian;
160}
161
162
163template<class Type>
166(
167 const areaScalarField& gamma,
169 const word& name
170)
171{
173 (
174 vf.mesh(),
175 vf.mesh().laplacianScheme(name)
176 ).ref().facLaplacian(gamma, vf);
177}
178
179template<class Type>
182(
183 const tmp<areaScalarField>& tgamma,
185 const word& name
186)
187{
189 (
190 fac::laplacian(tgamma(), vf, name)
191 );
192 tgamma.clear();
193 return Laplacian;
194}
195
196template<class Type>
199(
200 const areaScalarField& gamma,
202 const word& name
203)
204{
206 (
207 fac::laplacian(gamma, tvf(), name)
208 );
209 tvf.clear();
210 return Laplacian;
211}
212
213template<class Type>
216(
217 const tmp<areaScalarField>& tgamma,
219 const word& name
220)
221{
223 (
224 fac::laplacian(tgamma(), tvf(), name)
225 );
226 tgamma.clear();
227 tvf.clear();
228 return Laplacian;
229}
230
231
232template<class Type>
235(
236 const areaScalarField& gamma,
238)
239{
240 return fac::laplacian
241 (
242 gamma,
243 vf,
244 "laplacian(" + gamma.name() + ',' + vf.name() + ')'
245 );
246}
247
248template<class Type>
251(
252 const tmp<areaScalarField>& tgamma,
254)
255{
256 return fac::laplacian
257 (
258 tgamma,
259 vf,
260 "laplacian(" + tgamma().name() + ',' + vf.name() + ')'
261 );
262}
263
264template<class Type>
267(
268 const areaScalarField& gamma,
270)
271{
272 return fac::laplacian
273 (
274 gamma,
275 tvf,
276 "laplacian(" + gamma.name() + ',' + tvf().name() + ')'
277 );
278}
279
280template<class Type>
283(
284 const tmp<areaScalarField>& tgamma,
286)
287{
288 return fac::laplacian
289 (
290 tgamma,
291 tvf,
292 "laplacian(" + tgamma().name() + ',' + tvf().name() + ')'
293 );
294}
295
296
297template<class Type>
300(
301 const edgeScalarField& gamma,
303 const word& name
304)
305{
307 (
308 vf.mesh(),
309 vf.mesh().laplacianScheme(name)
310 ).ref().facLaplacian(gamma, vf);
311}
312
313template<class Type>
314tmp<GeometricField<Type, faPatchField, areaMesh>>
316(
317 const tmp<edgeScalarField>& tgamma,
318 const GeometricField<Type, faPatchField, areaMesh>& vf,
319 const word& name
320)
321{
322 tmp<GeometricField<Type, faPatchField, areaMesh>> Laplacian
323 (
324 fac::laplacian(tgamma(), vf, name)
325 );
326 tgamma.clear();
327 return Laplacian;
328}
329
330template<class Type>
331tmp<GeometricField<Type, faPatchField, areaMesh>>
333(
334 const edgeScalarField& gamma,
335 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tvf,
336 const word& name
337)
338{
339 tmp<GeometricField<Type, faPatchField, areaMesh>> Laplacian
340 (
341 fac::laplacian(gamma, tvf(), name)
342 );
343 tvf.clear();
344 return Laplacian;
345}
346
347template<class Type>
348tmp<GeometricField<Type, faPatchField, areaMesh>> laplacian
349(
350 const tmp<edgeScalarField>& tgamma,
351 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tvf,
352 const word& name
353)
354{
355 tmp<GeometricField<Type, faPatchField, areaMesh>> Laplacian
356 (
357 fac::laplacian(tgamma(), tvf(), name)
358 );
359 tgamma.clear();
360 tvf.clear();
361 return Laplacian;
362}
363
364
365template<class Type>
366tmp<GeometricField<Type, faPatchField, areaMesh>>
368(
369 const edgeScalarField& gamma,
370 const GeometricField<Type, faPatchField, areaMesh>& vf
371)
372{
373 return fac::laplacian
374 (
375 gamma,
376 vf,
377 "laplacian(" + gamma.name() + ',' + vf.name() + ')'
378 );
379}
380
381template<class Type>
382tmp<GeometricField<Type, faPatchField, areaMesh>>
384(
385 const tmp<edgeScalarField>& tgamma,
386 const GeometricField<Type, faPatchField, areaMesh>& vf
387)
388{
389 tmp<GeometricField<Type, faPatchField, areaMesh>> Laplacian
390 (
391 fac::laplacian(tgamma(), vf)
392 );
393 tgamma.clear();
394 return Laplacian;
395}
396
397template<class Type>
398tmp<GeometricField<Type, faPatchField, areaMesh>>
400(
401 const edgeScalarField& gamma,
402 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tvf
403)
404{
405 tmp<GeometricField<Type, faPatchField, areaMesh>> Laplacian
406 (
407 fac::laplacian(gamma, tvf())
408 );
409 tvf.clear();
410 return Laplacian;
411}
412
413template<class Type>
414tmp<GeometricField<Type, faPatchField, areaMesh>> laplacian
415(
416 const tmp<edgeScalarField>& tgamma,
417 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tvf
418)
419{
420 tmp<GeometricField<Type, faPatchField, areaMesh>> Laplacian
421 (
422 fac::laplacian(tgamma(), tvf())
423 );
424 tgamma.clear();
425 tvf.clear();
426 return Laplacian;
427}
428
429
430/*
431template<class Type>
432tmp<GeometricField<Type, faPatchField, areaMesh>>
433laplacian
434(
435 const areaTensorField& gamma,
436 const GeometricField<Type, faPatchField, areaMesh>& vf
437)
438{
439 const faMesh& mesh = tvf().mesh();
440 return fac::laplacian
441 (
442 (
443 mesh.Sf() & fac::interpolate(tgamma) & mesh.Sf()
444 )/sqr(mesh.magSf()),
445 tvf
446 );
447}
448
449template<class Type>
450tmp<GeometricField<Type, faPatchField, areaMesh>>
451laplacian
452(
453 const areaTensorField& gamma,
454 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tvf
455)
456{
457 tmp<GeometricField<Type, faPatchField, areaMesh>> Laplacian
458 (
459 fac::laplacian(gamma, tvf())
460 );
461 tvf.clear();
462 return Laplacian;
463}
464
465template<class Type>
466tmp<GeometricField<Type, faPatchField, areaMesh>>
467laplacian
468(
469 const tmp<areaTensorField>& tgamma,
470 const GeometricField<Type, faPatchField, areaMesh>& vf
471)
472{
473 tmp<GeometricField<Type, faPatchField, areaMesh>> Laplacian
474 (
475 fac::laplacian(tgamma(), vf)
476 );
477 tgamma.clear();
478 return Laplacian;
479}
480
481template<class Type>
482tmp<GeometricField<Type, faPatchField, areaMesh>>
483laplacian
484(
485 const tmp<areaTensorField>& tgamma,
486 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tvf
487)
488{
489 tmp<GeometricField<Type, faPatchField, areaMesh>> Laplacian
490 (
491 fac::laplacian(tgamma(), tvf())
492 );
493 tgamma.clear();
494 tvf.clear();
495 return Laplacian;
496}
497
498
499template<class Type>
500tmp<GeometricField<Type, faPatchField, areaMesh>> laplacian
501(
502 const edgeTensorField& gamma,
503 const GeometricField<Type, faPatchField, areaMesh>& vf
504)
505{
506 const faMesh& mesh = tvf().mesh();
507
508 return fac::laplacian
509 (
510 (
511 mesh.Sf() &
512 gamma
513 & mesh.Sf()
514 )/sqr(mesh.magSf()),
515 vf
516 );
517}
518
519template<class Type>
520tmp<GeometricField<Type, faPatchField, areaMesh>>
521laplacian
522(
523 const edgeTensorField& gamma,
524 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tvf
525)
526{
527 tmp<GeometricField<Type, faPatchField, areaMesh>> Laplacian
528 (
529 fac::laplacian(gamma, tvf())
530 );
531 tvf.clear();
532 return Laplacian;
533}
534
535template<class Type>
536tmp<GeometricField<Type, faPatchField, areaMesh>>
537laplacian
538(
539 const tmp<edgeTensorField>& tgamma,
540 const GeometricField<Type, faPatchField, areaMesh>& vf
541)
542{
543 tmp<GeometricField<Type, faPatchField, areaMesh>> Laplacian
544 (
545 fac::laplacian(tgamma(), tvf())
546 );
547 tgamma.clear();
548 return Laplacian;
549}
550
551template<class Type>
552tmp<GeometricField<Type, faPatchField, areaMesh>>
553laplacian
554(
555 const tmp<edgeTensorField>& tgamma,
556 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tvf
557)
558{
559 tmp<GeometricField<Type, faPatchField, areaMesh>> Laplacian
560 (
561 fac::laplacian(tgamma(), tvf())
562 );
563 tgamma.clear();
564 tvf.clear();
565 return Laplacian;
566}
567*/
568
569
570// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
571
572} // End namespace fac
573
574// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
575
576} // End namespace Foam
577
578// ************************************************************************* //
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A class for managing temporary objects.
Definition: tmp.H:65
void clear() const noexcept
Definition: tmpI.H:287
A class for handling words, derived from Foam::string.
Definition: word.H:68
const scalar gamma
Definition: EEqn.H:9
Calculate the laplacian of the given field.
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:47
Namespace for OpenFOAM.
GeometricField< scalar, faePatchField, edgeMesh > edgeScalarField
Definition: edgeFieldsFwd.H:63
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Calculate the second temporal derivative.