/// <summary>
        /// Constructor;
        /// </summary>
        /// <param name="Fields"></param>
        /// <param name="IterationResiduals"></param>
        /// <param name="LsTrk"></param>
        /// <param name="_ComputeOperatorMatrix">See <see cref="ComputeOperatorMatrix"/>.</param>
        /// <param name="_UpdateLevelset">See <see cref="UpdateLevelset"/>.</param>
        /// <param name="BDForder">
        /// The order of the BDF scheme from 1 to 6; in addition, 0 encodes Explicit Euler and -1 encodes Crank-Nicolson.
        /// </param>
        /// <param name="_LevelSetHandling"></param>
        /// <param name="_MassMatrixShapeandDependence"></param>
        /// <param name="_SpatialOperatorType"></param>
        /// <param name="_MassScale"></param>
        /// <param name="_AgglomerationThreshold"></param>
        /// <param name="_RKscheme"></param>
        /// <param name="useX">
        /// U dont want to know!
        /// </param>
        /// <param name="_MultigridSequence"></param>
        /// <param name="_CutCellQuadOrder">Order of quadrature in cut cells, required e.g. for <see cref="LevelSetTracker.GetXDGSpaceMetrics(SpeciesId[], int, int)"/></param>
        /// <param name="_SpId">Species to compute, actually a subset of <see cref="LevelSetTracker.SpeciesIdS"/></param>
        /// <param name="_MultigridOperatorConfig">
        /// Configuration of block-preconditioner, if null a default value is chosen.
        /// </param>
        public XdgRKTimestepping(DGField[] Fields,
                                 DGField[] IterationResiduals,
                                 LevelSetTracker LsTrk,
                                 DelComputeOperatorMatrix _ComputeOperatorMatrix,
                                 DelUpdateLevelset _UpdateLevelset,
                                 RungeKuttaScheme _RKscheme,
                                 LevelSetHandling _LevelSetHandling,
                                 MassMatrixShapeandDependence _MassMatrixShapeandDependence,
                                 SpatialOperatorType _SpatialOperatorType,
                                 IDictionary <SpeciesId, IEnumerable <double> > _MassScale,
                                 MultigridOperator.ChangeOfBasisConfig[][] _MultigridOperatorConfig,
                                 AggregationGridData[] _MultigridSequence,
                                 SpeciesId[] _SpId,
                                 int _CutCellQuadOrder,
                                 double _AgglomerationThreshold, bool useX,
                                 Control.NonLinearSolverConfig nonlinconfig,
                                 Control.LinearSolverConfig linearconfig) : base(nonlinconfig, linearconfig)
        {
            // check args, set internals
            // -------------------------

            if (Fields.Length != IterationResiduals.Length)
            {
                throw new ArgumentException("Expecting the same number of fields and residuals.");
            }
            for (int iFld = 0; iFld < Fields.Length; iFld++)
            {
                if (!Fields[iFld].Basis.Equals(IterationResiduals[iFld].Basis))
                {
                    throw new ArgumentException(string.Format("Mismatch between {0}-th basis of fields and residuals.", iFld));
                }
            }

            if (_MassScale != null)
            {
                if (!IEnumerableExtensions.SetEquals(_SpId, _MassScale.Keys))
                {
                    throw new ArgumentException();
                }
            }

            base.Residuals = new CoordinateVector(IterationResiduals);

            if (!(_RKscheme.IsExplicit || _RKscheme.IsDiagonallyImplicit))
            {
                throw new NotSupportedException("Only supporting explicit or diagonally implicit schemes.");
            }

            base.m_LsTrk = LsTrk;
            base.Config_LevelSetHandling             = _LevelSetHandling;
            base.Config_MassMatrixShapeandDependence = _MassMatrixShapeandDependence;
            base.Config_SpatialOperatorType          = _SpatialOperatorType;
            base.ComputeOperatorMatrix         = _ComputeOperatorMatrix;
            base.UpdateLevelset                = _UpdateLevelset;
            base.Config_MassScale              = _MassScale;
            base.Config_AgglomerationThreshold = _AgglomerationThreshold;
            this.m_RKscheme                    = _RKscheme.CloneAs();
            base.MultigridSequence             = _MultigridSequence;
            base.Config_SpeciesToCompute       = _SpId;
            base.Config_CutCellQuadratureOrder = _CutCellQuadOrder;
            if (_MultigridSequence == null || _MultigridSequence.Length < 1)
            {
                throw new ArgumentException("At least one grid level is required.");
            }

            m_CurrentState = new CoordinateVector(Fields);

            if (_MultigridOperatorConfig != null)
            {
                Config_MultigridOperator = _MultigridOperatorConfig;
            }
            else
            {
                SetConfig_MultigridOperator_Default(Fields);
            }

            base.CommonConfigurationChecks();

            // configure stack of level-set-tracker
            // ------------------------------------

            if (Config_LevelSetHandling == LevelSetHandling.None)
            {
                m_LsTrk.IncreaseHistoryLength(0);
            }
            else if (Config_LevelSetHandling == LevelSetHandling.LieSplitting ||
                     Config_LevelSetHandling == LevelSetHandling.StrangSplitting)
            {
                m_LsTrk.IncreaseHistoryLength(1);
            }
            else
            {
                m_LsTrk.IncreaseHistoryLength(m_RKscheme.Stages);
            }
            //m_LsTrk.IncreaseHistoryLength(1);

            // multigrid - init
            // ----------------

            InitMultigrid(Fields, useX);
        }
Exemple #2
0
        /*
         * /// <summary>
         * /// Legacy-Constructor for user-specified <see cref="DelComputeOperatorMatrix"/>
         * /// </summary>
         * public XdgTimestepping(
         *  DelComputeOperatorMatrix userComputeOperatorMatrix,
         *  IEnumerable<DGField> Fields,
         *  IEnumerable<DGField> IterationResiduals,
         *  TimeSteppingScheme __Scheme,
         *  DelUpdateLevelset _UpdateLevelset,
         *  LevelSetHandling _LevelSetHandling,
         *  MultigridOperator.ChangeOfBasisConfig[][] _MultigridOperatorConfig,
         *  AggregationGridData[] _MultigridSequence,
         *  double _AgglomerationThreshold,
         *  LinearSolverConfig LinearSolver, NonLinearSolverConfig NonLinearSolver) //
         * {
         *  this.Scheme = __Scheme;
         *  this.XdgOperator = op;
         *
         *  this.Parameters = op.InvokeParameterFactory(Fields);
         *
         *
         *  foreach (var f in Fields.Cat(IterationResiduals).Cat(Parameters)) {
         *      if (f != null && f is XDGField xf) {
         *          if (LsTrk == null) {
         *              LsTrk = xf.Basis.Tracker;
         *          } else {
         *              if (!object.ReferenceEquals(LsTrk, xf.Basis.Tracker))
         *                  throw new ArgumentException();
         *          }
         *      }
         *  }
         *  if (LsTrk == null)
         *      throw new ArgumentException("unable to get Level Set Tracker reference");
         *
         *  bool UseX = Fields.Any(f => f is XDGField) || IterationResiduals.Any(f => f is XDGField);
         *
         *  ConstructorCommon(op, UseX,
         *      Fields, this.Parameters, IterationResiduals,
         *      myDelComputeXOperatorMatrix,
         *      _UpdateLevelset,
         *      _LevelSetHandling,
         *      _MultigridOperatorConfig,
         *      _MultigridSequence,
         *      _AgglomerationThreshold,
         *      LinearSolver, NonLinearSolver);
         *
         * }
         */


        private void ConstructorCommon(
            ISpatialOperator op, bool UseX,
            IEnumerable <DGField> Fields, IEnumerable <DGField> __Parameters, IEnumerable <DGField> IterationResiduals,
            SpeciesId[] spcToCompute,
            DelUpdateLevelset _UpdateLevelset, LevelSetHandling _LevelSetHandling,
            MultigridOperator.ChangeOfBasisConfig[][] _MultigridOperatorConfig, AggregationGridData[] _MultigridSequence,
            double _AgglomerationThreshold,
            LinearSolverConfig LinearSolver, NonLinearSolverConfig NonLinearSolver) //
        {
            RungeKuttaScheme rksch;
            int bdfOrder;

            DecodeScheme(this.Scheme, out rksch, out bdfOrder);

            SpatialOperatorType _SpatialOperatorType = SpatialOperatorType.Nonlinear;


            int quadOrder = op.QuadOrderFunction(
                Fields.Select(f => f.Basis.Degree).ToArray(),
                Parameters.Select(f => f != null ? f.Basis.Degree : 0).ToArray(),
                IterationResiduals.Select(f => f.Basis.Degree).ToArray());

            // default solvers
            // ===============
            if (LinearSolver == null)
            {
                LinearSolver = new LinearSolverConfig()
                {
                    SolverCode = LinearSolverCode.automatic
                };
            }
            if (NonLinearSolver == null)
            {
                NonLinearSolver = new NonLinearSolverConfig()
                {
                    SolverCode = NonLinearSolverCode.Newton
                };
            }

            // default Multi-Grid
            // ==================

            if (_MultigridSequence == null)
            {
                _MultigridSequence = new[] { CoarseningAlgorithms.ZeroAggregation(this.GridDat) };
            }

            // default level-set treatment
            // ===========================

            if (_UpdateLevelset == null)
            {
                _UpdateLevelset = this.UpdateLevelsetWithNothing;
                if (_LevelSetHandling != LevelSetHandling.None)
                {
                    throw new ArgumentException($"If level-set handling is set to {_LevelSetHandling} (anything but {LevelSetHandling.None}) an updating routine must be specified.");
                }
            }

            // default multigrid operator config
            // =================================
            if (_MultigridOperatorConfig == null)
            {
                int NoOfVar = Fields.Count();
                _MultigridOperatorConfig    = new MultigridOperator.ChangeOfBasisConfig[0][];
                _MultigridOperatorConfig[0] = new MultigridOperator.ChangeOfBasisConfig[NoOfVar];
                for (int iVar = 0; iVar < NoOfVar; iVar++)
                {
                    _MultigridOperatorConfig[0][iVar] = new MultigridOperator.ChangeOfBasisConfig()
                    {
                        DegreeS  = new int[] { Fields.ElementAt(iVar).Basis.Degree },
                        mode     = MultigridOperator.Mode.Eye,
                        VarIndex = new int[] { iVar }
                    };
                }
            }

            // finally, create timestepper
            // ===========================

            if (bdfOrder > -1000)
            {
                m_BDF_Timestepper = new XdgBDFTimestepping(Fields, __Parameters, IterationResiduals,
                                                           LsTrk, true,
                                                           this.ComputeOperatorMatrix, op, _UpdateLevelset,
                                                           bdfOrder,
                                                           _LevelSetHandling,
                                                           MassMatrixShapeandDependence.IsTimeDependent,
                                                           _SpatialOperatorType,
                                                           _MultigridOperatorConfig, _MultigridSequence,
                                                           spcToCompute, quadOrder,
                                                           _AgglomerationThreshold, UseX,
                                                           NonLinearSolver,
                                                           LinearSolver);

                m_BDF_Timestepper.Config_AgglomerationThreshold = _AgglomerationThreshold;
            }
            else
            {
                m_RK_Timestepper = new XdgRKTimestepping(Fields.ToArray(), __Parameters, IterationResiduals.ToArray(),
                                                         LsTrk,
                                                         this.ComputeOperatorMatrix, op, _UpdateLevelset,
                                                         rksch,
                                                         _LevelSetHandling,
                                                         MassMatrixShapeandDependence.IsTimeDependent,
                                                         _SpatialOperatorType,
                                                         _MultigridOperatorConfig, _MultigridSequence,
                                                         spcToCompute, quadOrder,
                                                         _AgglomerationThreshold, UseX,
                                                         NonLinearSolver,
                                                         LinearSolver);

                m_RK_Timestepper.Config_AgglomerationThreshold = _AgglomerationThreshold;
            }
        }