LagrangianEqn.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "LagrangianEqn.H"
27 #include "GeometricField.H"
28 #include "toSubField.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class Type>
33 template<template<class> class PrimitiveField>
35 (
36  const word& name,
37  const tmp<LagrangianSubScalarField>& tDeltaT,
41 )
42 :
44  tDeltaT_(tDeltaT.valid() ? tDeltaT : tmp<LagrangianSubScalarField>()),
45  psiSubSub_(PsiRef<SubField>::ref(psi)),
46  psiSub_(PsiRef<Field>::ref(psi)),
47  psiPtr_(NullObjectPtr<LagrangianSubSubField<Type>>()),
48  deltaTSp0Ptr_(&deltaTSp0),
49  S0Ptr_(&S0),
50  deltaTSu(*this),
51  deltaTSp(*this),
52  Su(*this),
53  Sp(*this)
54 {}
55 
56 
57 template<class Type>
58 template<template<class> class PrimitiveField>
60 (
61  const word& name,
62  const LagrangianSubScalarField& deltaT,
66 )
67 :
69  (
70  name,
72  psi,
73  deltaTSp0,
74  S0
75  )
76 {}
77 
78 
79 template<class Type>
80 template<template<class> class PrimitiveField>
82 (
83  const word& name,
87 )
88 :
90 {}
91 
92 
93 template<class Type>
94 template<template<class> class PrimitiveField>
96 (
97  const tmp<LagrangianSubScalarField>& tDeltaT,
101 )
102 :
103  LagrangianEqnBase(word::null, psi.mesh()),
104  tDeltaT_(tDeltaT.valid() ? tDeltaT : tmp<LagrangianSubScalarField>()),
105  psiSubSub_(PsiRef<SubField>::ref(psi)),
106  psiSub_(PsiRef<Field>::ref(psi)),
107  psiPtr_(NullObjectPtr<LagrangianSubSubField<Type>>()),
108  deltaTSp0Ptr_(&deltaTSp0),
109  S0Ptr_(&S0),
110  deltaTSu(*this),
111  deltaTSp(*this),
112  Su(*this),
113  Sp(*this)
114 {}
115 
116 
117 template<class Type>
118 template<template<class> class PrimitiveField>
120 (
121  const LagrangianSubScalarField& deltaT,
125 )
126 :
127  LagrangianEqn(tmp<LagrangianSubScalarField>(deltaT), psi, deltaTSp0, S0)
128 {}
129 
130 
131 template<class Type>
132 template<template<class> class PrimitiveField>
134 (
138 )
139 :
140  LagrangianEqn(tmp<LagrangianSubScalarField>(), psi, deltaTSp0, S0)
141 {}
142 
143 
144 template<class Type>
146 (
147  const tmp<LagrangianSubScalarField>& tDeltaT,
151 )
152 :
153  LagrangianEqnBase(psi.name() + "Eqn", psi.mesh()),
154  tDeltaT_(tDeltaT.valid() ? tDeltaT : tmp<LagrangianSubScalarField>()),
155  psiSubSub_(psi),
156  psiSub_(NullObjectRef<LagrangianSubField<Type>>()),
157  psiPtr_(&psi),
158  deltaTSp0Ptr_(&deltaTSp0),
159  S0Ptr_(&S0),
160  deltaTSu(*this),
161  deltaTSp(*this),
162  Su(*this),
163  Sp(*this)
164 {}
165 
166 
167 template<class Type>
169 (
170  const LagrangianSubScalarField& deltaT,
174 )
175 :
176  LagrangianEqn(tmp<LagrangianSubScalarField>(deltaT), psi, deltaTSp0, S0)
177 {}
178 
179 
180 template<class Type>
182 (
186 )
187 :
188  LagrangianEqn(tmp<LagrangianSubScalarField>(), psi, deltaTSp0, S0)
189 {}
190 
191 
192 template<class Type>
194 :
195  tmp<LagrangianEqn<Type>>::refCount(),
196  LagrangianEqnBase(eqn),
197  tDeltaT_(eqn.tDeltaT_),
198  psiSubSub_(eqn.psiSubSub_),
199  psiSub_(eqn.psiSub_),
200  psiPtr_(eqn.psiPtr_),
201  deltaTSp0Ptr_(NullObjectPtr<LagrangianDynamicField<scalar>>()),
202  S0Ptr_(NullObjectPtr<LagrangianDynamicField<Type>>()),
203  deltaTSu(eqn.deltaTSu),
204  deltaTSp(eqn.deltaTSp),
205  Su(eqn.Su),
206  Sp(eqn.Sp)
207 {}
208 
209 
210 template<class Type>
212 :
213  tmp<LagrangianEqn<Type>>::refCount(),
214  LagrangianEqnBase(eqn),
215  tDeltaT_(move(eqn.tDeltaT_)),
216  psiSubSub_(eqn.psiSubSub_),
217  psiSub_(eqn.psiSub_),
218  psiPtr_(eqn.psiPtr_),
219  deltaTSp0Ptr_(eqn.deltaTSp0Ptr_),
220  S0Ptr_(eqn.S0Ptr_),
221  deltaTSu(move(eqn.deltaTSu)),
222  deltaTSp(move(eqn.deltaTSp)),
223  Su(move(eqn.Su)),
224  Sp(move(eqn.Sp))
225 {
226  eqn.deltaTSp0Ptr_ = NullObjectPtr<LagrangianDynamicField<scalar>>();
227  eqn.S0Ptr_ = NullObjectPtr<LagrangianDynamicField<Type>>();
228 }
229 
230 
231 template<class Type>
233 :
234  LagrangianEqnBase(tEqn()),
235  tDeltaT_(tEqn().tDeltaT_, tEqn.isTmp()),
236  psiSubSub_(tEqn().psiSubSub_),
237  psiSub_(tEqn().psiSub_),
238  psiPtr_(tEqn().psiPtr_),
239  deltaTSp0Ptr_
240  (
241  tEqn.isTmp()
242  ? tEqn().deltaTSp0Ptr_
244  ),
245  S0Ptr_
246  (
247  tEqn.isTmp()
248  ? tEqn().S0Ptr_
250  ),
251  deltaTSu
252  (
253  tEqn.isTmp()
254  ? LagrangianCoeff<Type, false>(tEqn.ref().deltaTSu, true)
255  : LagrangianCoeff<Type, false>(tEqn().deltaTSu)
256  ),
257  deltaTSp
258  (
259  tEqn.isTmp()
260  ? LagrangianCoeff<scalar, true>(tEqn.ref().deltaTSp, true)
261  : LagrangianCoeff<scalar, true>(tEqn().deltaTSp)
262  ),
263  Su
264  (
265  tEqn.isTmp()
266  ? LagrangianCoeff<Type, false>(tEqn.ref().Su, true)
267  : LagrangianCoeff<Type, false>(tEqn().Su)
268  ),
269  Sp
270  (
271  tEqn.isTmp()
272  ? LagrangianSp<Type>(tEqn.ref().Sp, true)
273  : LagrangianSp<Type>(tEqn().Sp)
274  )
275 {
276  if (tEqn.isTmp())
277  {
278  tEqn.ref().deltaTSp0Ptr_ =
279  NullObjectPtr<LagrangianDynamicField<scalar>>();
280  }
281  if (tEqn.isTmp())
282  {
283  tEqn.ref().S0Ptr_ =
284  NullObjectPtr<LagrangianDynamicField<Type>>();
285  }
286 
287  tEqn.clear();
288 }
289 
290 
291 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
292 
293 template<class Type>
295 {
296  if (notNull(deltaTSp0Ptr_))
297  {
299  (
300  mesh().sub(deltaTSp0Ptr_->internalField())
301  );
302 
303  deltaTSp0 = Zero;
304 
305  if (deltaTSp.valid()) deltaTSp0 = deltaTSp.S();
306  }
307 
308  if (notNull(S0Ptr_))
309  {
311  (
312  mesh().sub(S0Ptr_->internalField())
313  );
314 
315  mesh().group() == LagrangianGroup::none ? S0 = Zero : S0 *= -1;
316 
317  if (Su.valid()) S0 += Su.S();
318  if (Sp.valid()) S0 += Sp.Su(*this)->S();
319  }
320 }
321 
322 
323 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
324 
325 template<class Type>
327 (
328  const LagrangianEqn<Type>& eqn
329 )
330 {
331  return
333  (
334  notNull(eqn.psiSubSub_)
335  ? new LagrangianEqn<Type>(eqn.psiSubSub_)
336  : new LagrangianEqn<Type>(eqn.psiSub_)
337  );
338 }
339 
340 
341 template<class Type>
344 {
345  return
346  notNull(psiSubSub_)
347  ? tmp<LagrangianSubSubField<Type>>(psiSubSub_)
348  : toSubField(psiSub_);
349 }
350 
351 
352 template<class Type>
354 {
355  return
356  notNull(psiSubSub_)
357  ? static_cast<const regIOobject&>(psiSubSub_)
358  : static_cast<const regIOobject&>(psiSub_);
359 }
360 
361 
362 template<class Type>
364 {
365  return psiIo().name();
366 }
367 
368 
369 template<class Type>
370 template<template<class> class PrimitiveField>
372 (
374 ) const
375 {
376  return &static_cast<const regIOobject&>(psi) == &psiIo();
377 }
378 
379 
380 template<class Type>
382 {
383  if
384  (
385  tDeltaT_.valid()
386  && other.tDeltaT_.valid()
387  && &tDeltaT_() != &other.tDeltaT_()
388  )
389  {
391  << "Combining equations with different time-step fields"
392  << exit(FatalError);
393  }
394 
395  if (&psiIo() != &other.psiIo())
396  {
398  << "Combining equations with different fields"
399  << exit(FatalError);
400  }
401 
402  if (name_ == word::null || (isNull(psiPtr_) && notNull(other.psiPtr_)))
403  {
404  name_ = other.name_;
405  }
406  else if (name_ != other.name_)
407  {
408  name_ = word::null;
409  }
410 
411  if (!tDeltaT_.valid() && other.tDeltaT_.valid())
412  {
413  tDeltaT_ = tmp<LagrangianSubScalarField>(other.tDeltaT_);
414  }
415 
416  if (isNull(psiPtr_) && notNull(other.psiPtr_))
417  {
418  psiPtr_ = other.psiPtr_;
419  }
420 }
421 
422 
423 template<class Type>
425 (
426  const tmp<LagrangianEqn<Type>>& tOther
427 )
428 {
429  if (!tOther.isTmp()) return;
430 
431  if (notNull(tOther().deltaTSp0Ptr_))
432  {
433  if (notNull(deltaTSp0Ptr_) && deltaTSp0Ptr_ != tOther().deltaTSp0Ptr_)
434  {
436  << "Operating tmp equations with different previous "
437  << "implicit time coefficients" << exit(FatalError);
438  }
439 
440  if (isNull(deltaTSp0Ptr_)) deltaTSp0Ptr_ = tOther.ref().deltaTSp0Ptr_;
441 
442  tOther.ref().deltaTSp0Ptr_ =
443  NullObjectPtr<LagrangianDynamicField<scalar>>();
444  }
445 
446  if (notNull(tOther().S0Ptr_))
447  {
448  if (notNull(S0Ptr_) && S0Ptr_ != tOther().S0Ptr_)
449  {
451  << "Operating tmp equations with different previous "
452  << "sources" << exit(FatalError);
453  }
454 
455  if (isNull(S0Ptr_)) S0Ptr_ = tOther.ref().S0Ptr_;
456 
457  tOther.ref().S0Ptr_ =
458  NullObjectPtr<LagrangianDynamicField<Type>>();
459  }
460 
461  tOther.clear();
462 }
463 
464 
465 template<class Type>
468 {
469  if (!tDeltaT_.valid())
470  {
472  << "Cannot evaluate equation " << name_ << " for field "
473  << psiName() << " without a time-step" << exit(FatalError);
474  }
475 
477  (
479  );
480 
481  tResult.ref() += Su;
482  tResult.ref() *= tDeltaT_();
483  tResult.ref() += deltaTSu;
484 
485  return tResult;
486 }
487 
488 
489 template<class Type>
491 {
492  if (!tDeltaT_.valid())
493  {
495  << "Cannot evaluate equation " << name_ << " for field "
496  << psiName() << " without a time-step" << exit(FatalError);
497  }
498 
499  tmp<LagrangianSp<Type>> tResult(new LagrangianSp<Type>(*this));
500 
501  tResult.ref() += Sp;
502  tResult.ref() *= tDeltaT_();
503  tResult.ref() += deltaTSp;
504 
505  return tResult;
506 }
507 
508 
509 template<class Type>
512 {
513  if (!tDeltaT_.valid())
514  {
516  << "Cannot evaluate equation " << name_ << " for field "
517  << psiName() << " without a time-step" << exit(FatalError);
518  }
519 
521  (
523  );
524 
525  tResult.ref() += Su;
526  tResult.ref() += Sp.H();
527  tResult.ref() *= tDeltaT_();
528  tResult.ref() += deltaTSu;
529 
530  return tResult;
531 }
532 
533 
534 template<class Type>
537 {
538  if (!tDeltaT_.valid())
539  {
541  << "Cannot evaluate equation " << name_ << " for field "
542  << psiName() << " without a time-step" << exit(FatalError);
543  }
544 
546  (
548  );
549 
550  tResult.ref() += Sp.A();
551  tResult.ref() *= tDeltaT_();
552  tResult.ref() += deltaTSp;
553 
554  return tResult;
555 }
556 
557 
558 template<class Type>
559 void Foam::LagrangianEqn<Type>::solve(const bool final)
560 {
561  if (!tDeltaT_.valid())
562  {
564  << "Cannot solve equation " << name_ << " for field "
565  << psiName() << " without a time-step" << exit(FatalError);
566  }
567 
568  if (isNull(psiPtr_))
569  {
571  << "Cannot solve equation " << name_ << " for constant field "
572  << psiName() << exit(FatalError);
573  }
574 
575  if (!deltaTSp.valid() && !Sp.valid())
576  {
578  << "Cannot solve equation " << name_ << " for field "
579  << psiName() << " with no diagonal" << exit(FatalError);
580  }
581 
582  if (!deltaTSu.valid() && !Su.valid())
583  {
584  *psiPtr_ = Zero;
585  return;
586  }
587 
588  *psiPtr_ = - (allSu()()/allSp()());
589 
590  if (!final)
591  {
592  deltaTSp0Ptr_ = NullObjectPtr<LagrangianDynamicField<scalar>>();
593  S0Ptr_ = NullObjectPtr<LagrangianDynamicField<Type>>();
594  }
595 }
596 
597 
598 // * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
599 
600 template<class Type>
602 {
603  deltaTSu += other.deltaTSu;
604  deltaTSp += other.deltaTSp;
605  Su += other.Su;
606  Sp += other.Sp;
607 }
608 
609 
610 template<class Type>
612 (
613  const tmp<LagrangianEqn<Type>>& tOther
614 )
615 {
616  operator+=(tOther());
617  tOther.clear();
618 }
619 
620 
621 template<class Type>
623 {
624  deltaTSu -= other.deltaTSu;
625  deltaTSp -= other.deltaTSp;
626  Su -= other.Su;
627  Sp -= other.Sp;
628 }
629 
630 
631 template<class Type>
633 (
634  const tmp<LagrangianEqn<Type>>& tOther
635 )
636 {
637  operator-=(tOther());
638  tOther.clear();
639 }
640 
641 
642 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
643 
644 template<class Type>
646 (
647  const LagrangianEqn<Type>& eqn
648 )
649 {
650  tmp<LagrangianEqn<Type>> tResult(new LagrangianEqn<Type>(eqn));
651  tResult.ref().deltaTSu.negate();
652  tResult.ref().deltaTSp.negate();
653  tResult.ref().Su.negate();
654  tResult.ref().Sp.negate();
655  return tResult.ref();
656 }
657 
658 
659 template<class Type>
661 (
662  const tmp<LagrangianEqn<Type>>& tEqn
663 )
664 {
665  tmp<LagrangianEqn<Type>> tResult(-tEqn());
666  tEqn.clear();
667  return tResult;
668 }
669 
670 
671 #define LAGRANGIAN_EQN_EQN_OPERATOR(Op, EqOp) \
672  \
673  template<class Type> \
674  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
675  ( \
676  const LagrangianEqn<Type>& a, \
677  const LagrangianEqn<Type>& b \
678  ) \
679  { \
680  tmp<LagrangianEqn<Type>> tResult(new LagrangianEqn<Type>(a)); \
681  tResult.ref().op(b); \
682  tResult.ref() EqOp b; \
683  return tResult; \
684  }; \
685  \
686  template<class Type> \
687  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
688  ( \
689  const tmp<LagrangianEqn<Type>>& tA, \
690  const LagrangianEqn<Type>& b \
691  ) \
692  { \
693  tmp<LagrangianEqn<Type>> tResult(tA); \
694  tResult.ref().op(b); \
695  tResult.ref() EqOp b; \
696  return tResult; \
697  }; \
698  \
699  template<class Type> \
700  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
701  ( \
702  const tmp<LagrangianEqn<Type>>& tA, \
703  const tmp<LagrangianEqn<Type>>& tB \
704  ) \
705  { \
706  \
707  if (!tA.isTmp()) return tA() Op tB; \
708  if (!tB.isTmp()) return tA Op tB(); \
709  tmp<LagrangianEqn<Type>> tResult(tA Op tB()); \
710  tResult.ref().setPrevious(tB); \
711  return tResult; \
712  }
713 
714 #define LAGRANGIAN_COMMUTATIVE_EQN_EQN_OPERATOR(Op, EqOp) \
715  \
716  LAGRANGIAN_EQN_EQN_OPERATOR(Op, EqOp) \
717  \
718  template<class Type> \
719  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
720  ( \
721  const LagrangianEqn<Type>& a, \
722  const tmp<LagrangianEqn<Type>>& tB \
723  ) \
724  { \
725  return tB Op a; \
726  }
727 
728 #define LAGRANGIAN_NON_COMMUTATIVE_EQN_EQN_OPERATOR(Op, EqOp) \
729  \
730  LAGRANGIAN_EQN_EQN_OPERATOR(Op, EqOp) \
731  \
732  template<class Type> \
733  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
734  ( \
735  const LagrangianEqn<Type>& a, \
736  const tmp<LagrangianEqn<Type>>& tB \
737  ) \
738  { \
739  tmp<LagrangianEqn<Type>> tResult(a Op tB()); \
740  tResult.ref().setPrevious(tB); \
741  return tResult; \
742  }
743 
745 
747 
749 
750 #undef LAGRANGIAN_EQN_EQN_OPERATOR
751 #undef LAGRANGIAN_COMMUTATIVE_EQN_EQN_OPERATOR
752 #undef LAGRANGIAN_NON_COMMUTATIVE_EQN_EQN_OPERATOR
753 
754 
755 #define LAGRANGIAN_EQN_FIELD_OPERATOR(Op, EqOp, LagrangianSubField) \
756  \
757  template<class Type> \
758  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
759  ( \
760  const LagrangianEqn<Type>& eqn, \
761  const LagrangianSubField<Type>& field \
762  ) \
763  { \
764  tmp<LagrangianEqn<Type>> tResult(new LagrangianEqn<Type>(eqn)); \
765  tResult.ref().Su EqOp field; \
766  return tResult; \
767  } \
768  \
769  template<class Type> \
770  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
771  ( \
772  const tmp<LagrangianEqn<Type>>& tEqn, \
773  const LagrangianSubField<Type>& field \
774  ) \
775  { \
776  tmp<LagrangianEqn<Type>> tResult(tEqn); \
777  tResult.ref().Su EqOp field; \
778  return tResult; \
779  } \
780  \
781  template<class Type> \
782  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
783  ( \
784  const tmp<LagrangianEqn<Type>>& tEqn, \
785  const tmp<LagrangianSubField<Type>>& tField \
786  ) \
787  { \
788  tmp<LagrangianEqn<Type>> tResult(tEqn); \
789  tResult.ref().Su EqOp tField; \
790  return tResult; \
791  } \
792  \
793  template<class Type> \
794  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
795  ( \
796  const LagrangianEqn<Type>& eqn, \
797  const tmp<LagrangianSubField<Type>>& tField \
798  ) \
799  { \
800  tmp<LagrangianEqn<Type>> tResult(new LagrangianEqn<Type>(eqn)); \
801  tResult.ref().Su EqOp tField; \
802  return tResult; \
803  }
804 
805 #define LAGRANGIAN_FIELD_EQN_OPERATOR(Op, EqOp, LagrangianSubField) \
806  \
807  template<class Type> \
808  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
809  ( \
810  const LagrangianSubField<Type>& field, \
811  const LagrangianEqn<Type>& eqn \
812  ) \
813  { \
814  tmp<LagrangianEqn<Type>> tResult = \
815  LagrangianEqn<Type>::NewEmpty(eqn); \
816  tResult.ref() += field; \
817  tResult.ref() EqOp eqn; \
818  return tResult; \
819  } \
820  \
821  template<class Type> \
822  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
823  ( \
824  const LagrangianSubField<Type>& field, \
825  const tmp<LagrangianEqn<Type>>& tEqn \
826  ) \
827  { \
828  tmp<LagrangianEqn<Type>> tResult = \
829  LagrangianEqn<Type>::NewEmpty(tEqn()); \
830  tResult.ref() += field; \
831  tResult.ref() EqOp tEqn; \
832  return tResult; \
833  } \
834  \
835  template<class Type> \
836  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
837  ( \
838  const tmp<LagrangianSubField<Type>>& tField, \
839  const tmp<LagrangianEqn<Type>>& tEqn \
840  ) \
841  { \
842  tmp<LagrangianEqn<Type>> tResult = \
843  LagrangianEqn<Type>::NewEmpty(tEqn()); \
844  tResult.ref() += tField; \
845  tResult.ref() EqOp tEqn; \
846  return tResult; \
847  } \
848  \
849  template<class Type> \
850  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
851  ( \
852  const tmp<LagrangianSubField<Type>>& tField, \
853  const LagrangianEqn<Type>& eqn \
854  ) \
855  { \
856  tmp<LagrangianEqn<Type>> tResult = \
857  LagrangianEqn<Type>::NewEmpty(eqn); \
858  tResult.ref() += tField; \
859  tResult.ref() EqOp eqn; \
860  return tResult; \
861  }
862 
863 //- Addition
868 
869 //- Subtraction
874 
875 //- Set-equal-to
880 
881 #undef LAGRANGIAN_EQN_FIELD_OPERATOR
882 #undef LAGRANGIAN_FIELD_EQN_OPERATOR
883 
884 
885 #define LAGRANGIAN_EQN_SCALAR_FIELD_OPERATOR(Op, EqOp, LagrangianSubField) \
886  \
887  template<class Type> \
888  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
889  ( \
890  const LagrangianEqn<Type>& eqn, \
891  const LagrangianSubField<scalar>& field \
892  ) \
893  { \
894  tmp<LagrangianEqn<Type>> tResult(new LagrangianEqn<Type>(eqn)); \
895  tResult.ref().deltaTSu EqOp field; \
896  tResult.ref().deltaTSp EqOp field; \
897  tResult.ref().Su EqOp field; \
898  tResult.ref().Sp EqOp field; \
899  return tResult; \
900  } \
901  \
902  template<class Type> \
903  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
904  ( \
905  const tmp<LagrangianEqn<Type>>& tEqn, \
906  const LagrangianSubField<scalar>& field \
907  ) \
908  { \
909  tmp<LagrangianEqn<Type>> tResult(tEqn); \
910  tResult.ref().deltaTSu EqOp field; \
911  tResult.ref().deltaTSp EqOp field; \
912  tResult.ref().Su EqOp field; \
913  tResult.ref().Sp EqOp field; \
914  return tResult; \
915  } \
916  \
917  template<class Type> \
918  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
919  ( \
920  const tmp<LagrangianEqn<Type>>& tEqn, \
921  const tmp<LagrangianSubField<scalar>>& tField \
922  ) \
923  { \
924  tmp<LagrangianEqn<Type>> tResult(tEqn); \
925  tResult.ref().deltaTSu EqOp tField(); \
926  tResult.ref().deltaTSp EqOp tField(); \
927  tResult.ref().Su EqOp tField(); \
928  tResult.ref().Sp EqOp tField(); \
929  tField.clear(); \
930  return tResult; \
931  } \
932  \
933  template<class Type> \
934  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
935  ( \
936  const LagrangianEqn<Type>& eqn, \
937  const tmp<LagrangianSubField<scalar>>& tField \
938  ) \
939  { \
940  tmp<LagrangianEqn<Type>> tResult(new LagrangianEqn<Type>(eqn)); \
941  tResult.ref().deltaTSu EqOp tField(); \
942  tResult.ref().deltaTSp EqOp tField(); \
943  tResult.ref().Su EqOp tField(); \
944  tResult.ref().Sp EqOp tField(); \
945  tField.clear(); \
946  return tResult; \
947  }
948 
949 #define LAGRANGIAN_COMMUTATIVE_SCALAR_FIELD_EQN_OPERATOR(Op,LagrangianSubField)\
950  \
951  template<class Type> \
952  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
953  ( \
954  const LagrangianSubField<scalar>& field, \
955  const LagrangianEqn<Type>& eqn \
956  ) \
957  { \
958  return eqn Op field; \
959  } \
960  \
961  template<class Type> \
962  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
963  ( \
964  const LagrangianSubField<scalar>& field, \
965  const tmp<LagrangianEqn<Type>>& tEqn \
966  ) \
967  { \
968  return tEqn Op field; \
969  } \
970  \
971  template<class Type> \
972  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
973  ( \
974  const tmp<LagrangianSubField<scalar>>& tField, \
975  const tmp<LagrangianEqn<Type>>& tEqn \
976  ) \
977  { \
978  return tEqn Op tField; \
979  } \
980  \
981  template<class Type> \
982  Foam::tmp<Foam::LagrangianEqn<Type>> Foam::operator Op \
983  ( \
984  const tmp<LagrangianSubField<scalar>>& tField, \
985  const LagrangianEqn<Type>& eqn \
986  ) \
987  { \
988  return eqn Op tField; \
989  }
990 
991 //- Multiplication
996 
997 //- Division
1000 
1001 #undef LAGRANGIAN_EQN_SCALAR_FIELD_OPERATOR
1002 #undef LAGRANGIAN_COMMUTATIVE_SCALAR_FIELD_EQN_OPERATOR
1003 
1004 // ************************************************************************* //
#define LAGRANGIAN_COMMUTATIVE_SCALAR_FIELD_EQN_OPERATOR(Op, LagrangianSubField)
#define LAGRANGIAN_NON_COMMUTATIVE_EQN_EQN_OPERATOR(Op, EqOp)
#define LAGRANGIAN_COMMUTATIVE_EQN_EQN_OPERATOR(Op, EqOp)
#define LAGRANGIAN_EQN_SCALAR_FIELD_OPERATOR(Op, EqOp, LagrangianSubField)
#define LAGRANGIAN_EQN_FIELD_OPERATOR(Op, EqOp, LagrangianSubField)
#define LAGRANGIAN_FIELD_EQN_OPERATOR(Op, EqOp, LagrangianSubField)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Pre-declare SubField and related Field type.
Definition: Field.H:83
Generic GeometricField class.
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:131
Class to store a coefficient of a Lagrangian equation.
word name_
Symbolic name of the equation or equation term.
This class stores the coefficients of a Lagrangian equation, and facilitates solving that equation an...
Definition: LagrangianEqn.H:56
LagrangianCoeff< scalar, true > deltaTSp
Implicit time-coefficient.
Definition: LagrangianEqn.H:97
LagrangianSp< Type > Sp
Implicit coefficient.
void op(const LagrangianEqn< Type > &other)
Check the operation with another equation.
LagrangianCoeff< Type, false > Su
Explicit coefficient.
void solve(const bool final)
Solve.
tmp< LagrangianSp< Type > > allSp() const
Return the combined time and non-time implicit coefficient.
tmp< LagrangianCoeff< scalar, true > > allDiagonalSp() const
Return the combined time and non-time implicit diagonal.
const regIOobject & psiIo() const
Return the field IO.
~LagrangianEqn()
Destructor.
void operator+=(const LagrangianEqn< Type > &other)
Addition assignment.
tmp< LagrangianCoeff< Type, false > > allDiagonalSu() const
Return the combined time and non-time explicit diagonal.
tmp< LagrangianSubSubField< Type > > psi() const
Return the field.
static tmp< LagrangianEqn< Type > > NewEmpty(const LagrangianEqn< Type > &)
Construct from another equation, with empty coefficients.
LagrangianEqn(const word &name, const tmp< LagrangianSubScalarField > &tDeltaT, const LagrangianSubField< Type, PrimitiveField > &psi, LagrangianDynamicField< scalar > &deltaTSp0=NullObjectNonConstRef< LagrangianDynamicField< scalar >>(), LagrangianDynamicField< Type > &S0=NullObjectNonConstRef< LagrangianDynamicField< Type >>())
Construct for a const field and a tmp time-step with a name.
Definition: LagrangianEqn.C:35
void setPrevious(const tmp< LagrangianEqn< Type >> &tOther)
Set the previous field to that stored by another tmp equation.
LagrangianCoeff< Type, false > deltaTSu
Explicit time-coefficient.
Definition: LagrangianEqn.H:86
const word & psiName() const
Return the field name.
void operator-=(const LagrangianEqn< Type > &other)
Subtraction assignment.
bool isPsi(const LagrangianSubField< Type, PrimitiveField > &psi) const
Return whether the given field is that of the equation.
tmp< LagrangianCoeff< Type, false > > allSu() const
Return the combined time and non-time explicit coefficient.
Wrapper around LagrangianCoeff to specialise for the implicit coefficient. Trivial at present....
Definition: LagrangianSp.H:71
Pre-declare related SubField type.
Definition: SubField.H:63
Reference counter for various OpenFOAM components.
Definition: refCount.H:50
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:55
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
trAU ref().rename("rAU")
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const volScalarField & psi
bool valid(const PtrList< ModelType > &l)
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
DimensionedField< Type, LagrangianSubMesh, PrimitiveField > LagrangianSubField
DimensionedField< Type, LagrangianSubMesh, SubField > LagrangianSubSubField
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:64
T * NullObjectPtr()
Return pointer to the nullObject of type T.
Definition: nullObjectI.H:45
const T & NullObjectRef()
Return const reference to the nullObject of type T.
Definition: nullObjectI.H:27
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:58
error FatalError
void operator+=(fvMatrix< Type > &fvEqn, const CarrierEqn< Type > &cEqn)
Add to a finite-volume equation.
Definition: CarrierEqn.C:71
void operator-=(fvMatrix< Type > &fvEqn, const CarrierEqn< Type > &cEqn)
Subtract from a finite-volume equation.
Definition: CarrierEqn.C:94
tmp< DimensionedField< Type, GeoMesh, SubField > > toSubField(const DimensionedField< Type, GeoMesh, Field > &)
Return a temporary sub-field from a reference to a field.
Functions to cast/convert dimensioned field references and temporaries based on a primitive field to ...