/// <summary> /// Difference between <paramref name="f1"/> and <paramref name="f0"/> /// in the L2Norm. /// </summary> /// <param name="f1"></param> /// <param name="f0"></param> /// <param name="ResidualKey"> /// Optional parameter for key of residual. /// If null, the name of <paramref name="f1"/> will be used as key. /// </param> public void ComputeDifference(DGField f1, DGField f0, string ResidualKey = null) { if (f1.Basis.Degree != f0.Basis.Degree) { throw new ArgumentException("Fields must have equal degrees!"); } string key = GetDictionaryKey(Residualtype.L2Norm, f1.Identification, ResidualKey); SinglePhaseField diff = new SinglePhaseField(f1.Basis); diff.Acc(1.0, f1); diff.Acc(-1.0, f0); double res = diff.L2Norm(); m_Residuals[key] = res; }
protected override double RunSolverOneStep(int TimestepNo, double phystime, double dt) { Console.WriteLine(" Timestep # " + TimestepNo + ", phystime = " + phystime); //phystime = 1.8; LsUpdate(phystime); // operator-matrix assemblieren MsrMatrix OperatorMatrix = new MsrMatrix(u.Mapping, u.Mapping); double[] Affine = new double[OperatorMatrix.RowPartitioning.LocalLength]; MultiphaseCellAgglomerator Agg; MassMatrixFactory Mfact; // Agglomerator setup int quadOrder = Op.QuadOrderFunction(new int[] { u.Basis.Degree }, new int[0], new int[] { u.Basis.Degree }); //Agg = new MultiphaseCellAgglomerator(new CutCellMetrics(MomentFittingVariant, quadOrder, LsTrk, ), this.THRESHOLD, false); Agg = LsTrk.GetAgglomerator(new SpeciesId[] { LsTrk.GetSpeciesId("B") }, quadOrder, this.THRESHOLD); Console.WriteLine("Inter-Process agglomeration? " + Agg.GetAgglomerator(LsTrk.GetSpeciesId("B")).AggInfo.InterProcessAgglomeration); if (this.THRESHOLD > 0.01) { TestAgglomeration_Extraploation(Agg); TestAgglomeration_Projection(quadOrder, Agg); } // operator matrix assembly Op.ComputeMatrixEx(LsTrk, u.Mapping, null, u.Mapping, OperatorMatrix, Affine, false, 0.0, true, Agg.CellLengthScales, LsTrk.GetSpeciesId("B")); Agg.ManipulateMatrixAndRHS(OperatorMatrix, Affine, u.Mapping, u.Mapping); // mass matrix factory Mfact = LsTrk.GetXDGSpaceMetrics(new SpeciesId[] { LsTrk.GetSpeciesId("B") }, quadOrder, 1).MassMatrixFactory;// new MassMatrixFactory(u.Basis, Agg); // Mass matrix/Inverse Mass matrix //var MassInv = Mfact.GetMassMatrix(u.Mapping, new double[] { 1.0 }, true, LsTrk.GetSpeciesId("B")); var Mass = Mfact.GetMassMatrix(u.Mapping, new double[] { 1.0 }, false, LsTrk.GetSpeciesId("B")); Agg.ManipulateMatrixAndRHS(Mass, default(double[]), u.Mapping, u.Mapping); var MassInv = Mass.InvertBlocks(OnlyDiagonal: true, Subblocks: true, ignoreEmptyBlocks: true, SymmetricalInversion: false); // test that operator depends only on B-species values double DepTest = LsTrk.Regions.GetSpeciesSubGrid("B").TestMatrixDependency(OperatorMatrix, u.Mapping, u.Mapping); Console.WriteLine("Matrix dependency test: " + DepTest); Assert.LessOrEqual(DepTest, 0.0); // diagnostic output Console.WriteLine("Number of Agglomerations (all species): " + Agg.TotalNumberOfAgglomerations); Console.WriteLine("Number of Agglomerations (species 'B'): " + Agg.GetAgglomerator(LsTrk.GetSpeciesId("B")).AggInfo.SourceCells.NoOfItemsLocally.MPISum()); // operator auswerten: double[] x = new double[Affine.Length]; BLAS.daxpy(x.Length, 1.0, Affine, 1, x, 1); OperatorMatrix.SpMVpara(1.0, u.CoordinateVector, 1.0, x); MassInv.SpMV(1.0, x, 0.0, du_dx.CoordinateVector); Agg.GetAgglomerator(LsTrk.GetSpeciesId("B")).Extrapolate(du_dx.Mapping); // markieren, wo ueberhaupt A und B sind Bmarker.AccConstant(1.0, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); Amarker.AccConstant(+1.0, LsTrk.Regions.GetSpeciesSubGrid("A").VolumeMask); Xmarker.AccConstant(+1.0, LsTrk.Regions.GetSpeciesSubGrid("X").VolumeMask); // compute error ERR.Clear(); ERR.Acc(1.0, du_dx_Exact, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); ERR.Acc(-1.0, du_dx, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); double L2Err = ERR.L2Norm(LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); Console.WriteLine("L2 Error: " + L2Err); XERR.Clear(); XERR.GetSpeciesShadowField("B").Acc(1.0, ERR, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); double xL2Err = XERR.L2Norm(); Console.WriteLine("L2 Error (in XDG space): " + xL2Err); // check error if (this.THRESHOLD > 0.01) { // without agglomeration, the error in very tiny cut-cells may be large over the whole cell // However, the error in the XDG-space should be small under all circumstances Assert.LessOrEqual(L2Err, 1.0e-6); } Assert.LessOrEqual(xL2Err, 1.0e-6); bool IsPassed = ((L2Err <= 1.0e-6 || this.THRESHOLD <= 0.01) && xL2Err <= 1.0e-7); if (IsPassed) { Console.WriteLine("Test PASSED"); } else { Console.WriteLine("Test FAILED: check errors."); } // return/Ende base.NoOfTimesteps = 17; //base.NoOfTimesteps = 2; dt = 0.3; return(dt); }
/// <summary> /// One Iteration of the ReInitialization /// Operators must be built first /// </summary> /// <param name="ChangeRate"> /// L2-Norm of the Change-Rate in the level set in this reinit step /// </param> /// <param name="Restriction"> /// The subgrid, on which the ReInit is performed /// </param> /// <param name="IncludingInterface"> /// !! Not yet functional !! /// True, if the subgrid contains the interface, this causes all external edges of the subgrid to be treated as boundaries /// False, for the rest of the domain, thus the flux to the adjacent cells wil be evaluated /// </param> /// <returns></returns> public void ReInitSingleStep(out double ChangeRate, SubGrid Restriction = null, bool IncludingInterface = true) { if (!IncludingInterface) { throw new NotImplementedException("Untested, not yet functional!"); } using (new FuncTrace()) { /// Init Residuals Residual.Clear(); Residual.Acc(1.0, Phi); OldPhi.Clear(); OldPhi.Acc(1.0, Phi); NewPhi.Clear(); NewPhi.Acc(1.0, Phi); CellMask RestrictionMask = Restriction == null ? null : Restriction.VolumeMask; //if (Control.Upwinding && UpdateDirection && IterationCounter % 10 == 0) { if (false && Control.Upwinding && UpdateDirection) { //if (Control.Upwinding && UpdateDirection) { UpdateBulkMatrix(Restriction); } UpdateDirection = false; // RHS part RHSField.CoordinateVector.Clear(); //Operator_RHS.Evaluate(NewPhi.Mapping, RHSField.Mapping); Operator_RHS.Evaluate(double.NaN, IncludingInterface ? Restriction : null, IncludingInterface ? SubGridBoundaryModes.BoundaryEdge : SubGridBoundaryModes.InnerEdge, ArrayTools.Cat(new DGField[] { Phi }, parameterFields, new DGField[] { RHSField })); #if DEBUG RHSField.CheckForNanOrInf(); #endif // solve // ===== double[] RHS = OpAffine.CloneAs(); RHS.ScaleV(-1.0); RHS.AccV(1.0, RHSField.CoordinateVector); SolverResult Result; if (Restriction != null) { SubRHS.Clear(); SubSolution.Clear(); SubRHS.AccV(1.0, RHS, default(int[]), SubVecIdx); SubSolution.AccV(1.0, NewPhi.CoordinateVector, default(int[]), SubVecIdx); Result = slv.Solve(SubSolution, SubRHS); NewPhi.Clear(RestrictionMask); NewPhi.CoordinateVector.AccV(1.0, SubSolution, SubVecIdx, default(int[])); } else { Result = slv.Solve(NewPhi.CoordinateVector, RHS); } #if Debug OpMatrix.SpMV(-1.0, NewPhi.CoordinateVector, 1.0, RHS); Console.WriteLine("LinearSolver: {0} Iterations, Converged={1}, Residual = {2} ", Result.NoOfIterations, Result.Converged, RHS.L2Norm()); #endif // Apply underrelaxation Phi.Clear(RestrictionMask); Phi.Acc(1 - underrelaxation, OldPhi, RestrictionMask); Phi.Acc(underrelaxation, NewPhi, RestrictionMask); Residual.Acc(-1.0, Phi, RestrictionMask); ChangeRate = Residual.L2Norm(RestrictionMask); //Calculate LevelSetGradient.Clear(); LevelSetGradient.Gradient(1.0, Phi, RestrictionMask); //LevelSetGradient.GradientByFlux(1.0, Phi); MeanLevelSetGradient.Clear(); MeanLevelSetGradient.AccLaidBack(1.0, LevelSetGradient, RestrictionMask); if (Control.Upwinding) { //RestrictionMask.GetBitMask(); for (int i = 0; i < MeanLevelSetGradient.CoordinateVector.Length; i++) { NewDirection[i] = Math.Sign(MeanLevelSetGradient.CoordinateVector[i]); //NewDirection[i] = MeanLevelSetGradient.CoordinateVector[i]; OldDirection[i] -= NewDirection[i]; } double MaxDiff = OldDirection.L2Norm(); //if (MaxDiff > 1E-20 && IterationCounter % 10 == 0 ) { //if (MaxDiff > 1E-20) { // Console.WriteLine("Direction Values differ by {0}", MaxDiff); if (MaxDiff > 0.2) { //UpdateDirection = true; //Console.WriteLine("Direction Values differ by {0} => Updating ReInit-Matrix", MaxDiff); } ; //} //Console.WriteLine("HACK!!! Updating Upwind Matrix everytime!"); //UpdateDirection = true; // Reset Value OldDirection.Clear(); OldDirection.AccV(1.0, NewDirection); } } }
/// <summary> /// Computes the new level set field at time <paramref name="Phystime"/> + <paramref name="dt"/>. /// This is a 'driver function' which provides a universal interface to the various level set evolution algorithms. /// It also acts as a callback to the time stepper (see <see cref="m_BDF_Timestepper"/> resp. <see cref="m_RK_Timestepper"/>), /// i.e. it matches the signature of /// <see cref="BoSSS.Solution.XdgTimestepping.DelUpdateLevelset"/>. /// </summary> /// <param name="Phystime"></param> /// <param name="dt"></param> /// <param name="CurrentState"> /// The current solution (velocity and pressure), since the complete solution is provided by the time stepper, /// only the velocity components(supposed to be at the beginning) are used. /// </param> /// <param name="underrelax"> /// </param> public double DelUpdateLevelSet(DGField[] CurrentState, double Phystime, double dt, double underrelax, bool incremental) { using (new FuncTrace()) { //dt *= underrelax; int D = base.Grid.SpatialDimension; int iTimestep = hack_TimestepIndex; DGField[] EvoVelocity = CurrentState.GetSubVector(0, D); // ======================================================== // Backup old level-set, in order to compute the residual // ======================================================== SinglePhaseField LsBkUp = new SinglePhaseField(this.LevSet.Basis); LsBkUp.Acc(1.0, this.LevSet); CellMask oldCC = LsTrk.Regions.GetCutCellMask(); // ==================================================== // set evolution velocity, but only on the CUT-cells // ==================================================== #region Calculate density averaged Velocity for each cell ConventionalDGField[] meanVelocity = GetMeanVelocityFromXDGField(EvoVelocity); #endregion // =================================================================== // backup interface properties (mass conservation, surface changerate) // =================================================================== #region backup interface props double oldSurfVolume = 0.0; double oldSurfLength = 0.0; double SurfChangerate = 0.0; if (this.Control.CheckInterfaceProps) { oldSurfVolume = XNSEUtils.GetSpeciesArea(this.LsTrk, LsTrk.GetSpeciesId("A")); oldSurfLength = XNSEUtils.GetInterfaceLength(this.LsTrk); SurfChangerate = EnergyUtils.GetSurfaceChangerate(this.LsTrk, meanVelocity, this.m_HMForder); } #endregion // ==================================================== // perform level-set evolution // ==================================================== #region level-set evolution // set up for Strang splitting SinglePhaseField DGLevSet_old; if (incremental) { DGLevSet_old = this.DGLevSet.Current.CloneAs(); } else { DGLevSet_old = this.DGLevSet[0].CloneAs(); } // set up for underrelaxation SinglePhaseField DGLevSet_oldIter = this.DGLevSet.Current.CloneAs(); //PlotCurrentState(hack_Phystime, new TimestepNumber(new int[] { hack_TimestepIndex, 0 }), 2); // actual evolution switch (this.Control.Option_LevelSetEvolution) { case LevelSetEvolution.None: throw new ArgumentException("illegal call"); case LevelSetEvolution.FastMarching: { NarrowMarchingBand.Evolve_Mk2( dt, this.LsTrk, DGLevSet_old, this.DGLevSet.Current, this.DGLevSetGradient, meanVelocity, this.ExtensionVelocity.Current.ToArray(), //new DGField[] { LevSetSrc }, this.m_HMForder, iTimestep); //FastMarchReinitSolver = new FastMarchReinit(DGLevSet.Current.Basis); //CellMask Accepted = LsTrk.Regions.GetCutCellMask(); //CellMask ActiveField = LsTrk.Regions.GetNearFieldMask(1); //CellMask NegativeField = LsTrk.Regions.GetSpeciesMask("A"); //FastMarchReinitSolver.FirstOrderReinit(DGLevSet.Current, Accepted, NegativeField, ActiveField); break; } case LevelSetEvolution.Fourier: { Fourier_Timestepper.moveLevelSet(dt, meanVelocity); if (incremental) { Fourier_Timestepper.updateFourierLevSet(); } Fourier_LevSet.ProjectToDGLevelSet(this.DGLevSet.Current, this.LsTrk); break; } case LevelSetEvolution.Prescribed: { this.DGLevSet.Current.Clear(); this.DGLevSet.Current.ProjectField(1.0, this.Control.Phi.Vectorize(Phystime + dt)); break; } case LevelSetEvolution.ScalarConvection: { var LSM = new LevelSetMover(EvoVelocity, this.ExtensionVelocity, this.LsTrk, XVelocityProjection.CutCellVelocityProjectiontype.L2_plain, this.DGLevSet, this.BcMap); int check1 = this.ExtensionVelocity.PushCount; int check2 = this.DGLevSet.PushCount; this.DGLevSet[1].Clear(); this.DGLevSet[1].Acc(1.0, DGLevSet_old); LSM.Advect(dt); if (check1 != this.ExtensionVelocity.PushCount) { throw new ApplicationException(); } if (check2 != this.DGLevSet.PushCount) { throw new ApplicationException(); } break; } case LevelSetEvolution.ExtensionVelocity: { DGLevSetGradient.Clear(); DGLevSetGradient.Gradient(1.0, DGLevSet.Current); ExtVelMover.Advect(dt); // Fast Marching: Specify the Domains first // Perform Fast Marching only on the Far Field if (this.Control.AdaptiveMeshRefinement) { int NoCells = ((GridData)this.GridData).Cells.Count; BitArray Refined = new BitArray(NoCells); for (int j = 0; j < NoCells; j++) { if (((GridData)this.GridData).Cells.GetCell(j).RefinementLevel > 0) { Refined[j] = true; } } CellMask Accepted = new CellMask(this.GridData, Refined); CellMask AcceptedNeigh = Accepted.AllNeighbourCells(); Accepted = Accepted.Union(AcceptedNeigh); CellMask ActiveField = Accepted.Complement(); CellMask NegativeField = LsTrk.Regions.GetSpeciesMask("A"); FastMarchReinitSolver.FirstOrderReinit(DGLevSet.Current, Accepted, NegativeField, ActiveField); } else { CellMask Accepted = LsTrk.Regions.GetNearFieldMask(1); CellMask ActiveField = Accepted.Complement(); CellMask NegativeField = LsTrk.Regions.GetSpeciesMask("A"); FastMarchReinitSolver.FirstOrderReinit(DGLevSet.Current, Accepted, NegativeField, ActiveField); } //SubGrid AcceptedGrid = new SubGrid(Accepted); //ReInitPDE.ReInitialize(Restriction: AcceptedGrid); //CellMask ActiveField = Accepted.Complement(); //CellMask NegativeField = LsTrk.Regions.GetSpeciesMask("A"); //FastMarchReinitSolver.FirstOrderReinit(DGLevSet.Current, Accepted, NegativeField, ActiveField); //ReInitPDE.ReInitialize(); break; } default: throw new ApplicationException(); } // performing underrelaxation if (underrelax < 1.0) { this.DGLevSet.Current.Scale(underrelax); this.DGLevSet.Current.Acc((1.0 - underrelax), DGLevSet_oldIter); } //PlotCurrentState(hack_Phystime, new TimestepNumber(new int[] { hack_TimestepIndex, 1 }), 2); #endregion // ====================== // postprocessing // ======================= if (this.Control.ReInitPeriod > 0 && hack_TimestepIndex % this.Control.ReInitPeriod == 0) { Console.WriteLine("Filtering DG-LevSet"); SinglePhaseField FiltLevSet = new SinglePhaseField(DGLevSet.Current.Basis); FiltLevSet.AccLaidBack(1.0, DGLevSet.Current); Filter(FiltLevSet, 2, oldCC); DGLevSet.Current.Clear(); DGLevSet.Current.Acc(1.0, FiltLevSet); Console.WriteLine("FastMarchReInit performing FirstOrderReInit"); FastMarchReinitSolver = new FastMarchReinit(DGLevSet.Current.Basis); CellMask Accepted = LsTrk.Regions.GetCutCellMask(); CellMask ActiveField = LsTrk.Regions.GetNearFieldMask(1); CellMask NegativeField = LsTrk.Regions.GetSpeciesMask("A"); FastMarchReinitSolver.FirstOrderReinit(DGLevSet.Current, Accepted, NegativeField, ActiveField); } #region ensure continuity // make level set continuous CellMask CC = LsTrk.Regions.GetCutCellMask4LevSet(0); CellMask Near1 = LsTrk.Regions.GetNearMask4LevSet(0, 1); CellMask PosFF = LsTrk.Regions.GetLevelSetWing(0, +1).VolumeMask; ContinuityEnforcer.MakeContinuous(this.DGLevSet.Current, this.LevSet, Near1, PosFF); if (this.Control.Option_LevelSetEvolution == LevelSetEvolution.FastMarching) { CellMask Nearband = Near1.Union(CC); this.DGLevSet.Current.Clear(Nearband); this.DGLevSet.Current.AccLaidBack(1.0, this.LevSet, Nearband); //ContinuityEnforcer.SetFarField(this.DGLevSet.Current, Near1, PosFF); } //PlotCurrentState(hack_Phystime, new TimestepNumber(new int[] { hack_TimestepIndex, 2 }), 2); #endregion for (int d = 0; d < D; d++) { this.XDGvelocity.Velocity[d].UpdateBehaviour = BehaveUnder_LevSetMoovement.AutoExtrapolate; } // =============== // tracker update // =============== this.LsTrk.UpdateTracker(Phystime + dt, incremental: true); // update near field (in case of adaptive mesh refinement) if (this.Control.AdaptiveMeshRefinement && this.Control.Option_LevelSetEvolution == LevelSetEvolution.FastMarching) { Near1 = LsTrk.Regions.GetNearMask4LevSet(0, 1); PosFF = LsTrk.Regions.GetLevelSetWing(0, +1).VolumeMask; ContinuityEnforcer.SetFarField(this.DGLevSet.Current, Near1, PosFF); ContinuityEnforcer.SetFarField(this.LevSet, Near1, PosFF); } // ================================================================== // check interface properties (mass conservation, surface changerate) // ================================================================== if (this.Control.CheckInterfaceProps) { double currentSurfVolume = XNSEUtils.GetSpeciesArea(this.LsTrk, LsTrk.GetSpeciesId("A")); double massChange = ((currentSurfVolume - oldSurfVolume) / oldSurfVolume) * 100; Console.WriteLine("Change of mass = {0}%", massChange); double currentSurfLength = XNSEUtils.GetInterfaceLength(this.LsTrk); double actualSurfChangerate = (currentSurfLength - oldSurfLength) / dt; Console.WriteLine("Interface divergence = {0}", SurfChangerate); Console.WriteLine("actual surface changerate = {0}", actualSurfChangerate); } // ================== // compute residual // ================== var newCC = LsTrk.Regions.GetCutCellMask(); LsBkUp.Acc(-1.0, this.LevSet); double LevSetResidual = LsBkUp.L2Norm(newCC.Union(oldCC)); return(LevSetResidual); } }
/// <summary> /// Single run of the solver /// </summary> protected override double RunSolverOneStep(int TimestepNo, double phystime, double dt) { using (new FuncTrace()) { if (Control.ExactSolution_provided) { Tex.Clear(); Tex.ProjectField(this.Control.InitialValues_Evaluators["Tex"]); RHS.Clear(); RHS.ProjectField(this.Control.InitialValues_Evaluators["RHS"]); } if (Control.AdaptiveMeshRefinement == false) { base.NoOfTimesteps = -1; if (TimestepNo > 1) { throw new ApplicationException("steady-state-equation."); } base.TerminationKey = true; } // Update matrices // --------------- UpdateMatrices(); // call solver // ----------- double mintime, maxtime; bool converged; int NoOfIterations; ClassicSolve(out mintime, out maxtime, out converged, out NoOfIterations); Console.WriteLine("finished; " + NoOfIterations + " iterations."); Console.WriteLine("converged? " + converged); Console.WriteLine("Timespan: " + mintime + " to " + maxtime + " seconds"); base.QueryHandler.ValueQuery("minSolRunT", mintime, true); base.QueryHandler.ValueQuery("maxSolRunT", maxtime, true); base.QueryHandler.ValueQuery("Conv", converged ? 1.0 : 0.0, true); base.QueryHandler.ValueQuery("NoIter", NoOfIterations, true); base.QueryHandler.ValueQuery("NoOfCells", this.GridData.CellPartitioning.TotalLength, true); base.QueryHandler.ValueQuery("DOFs", T.Mapping.TotalLength, true); base.QueryHandler.ValueQuery("BlockSize", T.Basis.Length, true); if (base.Control.ExactSolution_provided) { Error.Clear(); Error.AccLaidBack(1.0, Tex); Error.AccLaidBack(-1.0, T); double L2_ERR = Error.L2Norm(); Console.WriteLine("\t\tL2 error on " + this.Grid.NumberOfCells + ": " + L2_ERR); base.QueryHandler.ValueQuery("SolL2err", L2_ERR, true); } // evaluate residual // ================= { this.ResiualKP1.Clear(); if (this.Control.InitialValues_Evaluators.ContainsKey("RHS")) { this.ResiualKP1.ProjectField(this.Control.InitialValues_Evaluators["RHS"]); } var ev = this.LapaceIp.GetEvaluator(T.Mapping, ResiualKP1.Mapping); ev.Evaluate(-1.0, 1.0, ResiualKP1.CoordinateVector); } // return // ====== return(0.0); } }
protected override double RunSolverOneStep(int TimestepNo, double phystime, double dt) { Console.WriteLine(" Timestep # " + TimestepNo + ", phystime = " + phystime); //phystime = 1.8; LsUpdate(phystime); // operator-matrix assemblieren MsrMatrix OperatorMatrix = new MsrMatrix(u.Mapping, u.Mapping); double[] Affine = new double[OperatorMatrix.RowPartitioning.LocalLength]; // Agglomerator setup MultiphaseCellAgglomerator Agg = LsTrk.GetAgglomerator(new SpeciesId[] { LsTrk.GetSpeciesId("B") }, QuadOrder, this.THRESHOLD); // plausibility of cell length scales if (SER_PAR_COMPARISON) { TestLengthScales(QuadOrder, TimestepNo); } Console.WriteLine("Inter-Process agglomeration? " + Agg.GetAgglomerator(LsTrk.GetSpeciesId("B")).AggInfo.InterProcessAgglomeration); if (this.THRESHOLD > 0.01) { TestAgglomeration_Extraploation(Agg); TestAgglomeration_Projection(QuadOrder, Agg); } CheckExchange(true); CheckExchange(false); // operator matrix assembly XSpatialOperatorMk2.XEvaluatorLinear mtxBuilder = Op.GetMatrixBuilder(base.LsTrk, u.Mapping, null, u.Mapping); mtxBuilder.time = 0.0; mtxBuilder.ComputeMatrix(OperatorMatrix, Affine); Agg.ManipulateMatrixAndRHS(OperatorMatrix, Affine, u.Mapping, u.Mapping); // mass matrix factory var Mfact = LsTrk.GetXDGSpaceMetrics(new SpeciesId[] { LsTrk.GetSpeciesId("B") }, QuadOrder, 1).MassMatrixFactory;// new MassMatrixFactory(u.Basis, Agg); // Mass matrix/Inverse Mass matrix //var MassInv = Mfact.GetMassMatrix(u.Mapping, new double[] { 1.0 }, true, LsTrk.GetSpeciesId("B")); var Mass = Mfact.GetMassMatrix(u.Mapping, new double[] { 1.0 }, false, LsTrk.GetSpeciesId("B")); Agg.ManipulateMatrixAndRHS(Mass, default(double[]), u.Mapping, u.Mapping); var MassInv = Mass.InvertBlocks(OnlyDiagonal: true, Subblocks: true, ignoreEmptyBlocks: true, SymmetricalInversion: false); // test that operator depends only on B-species values double DepTest = LsTrk.Regions.GetSpeciesSubGrid("B").TestMatrixDependency(OperatorMatrix, u.Mapping, u.Mapping); Console.WriteLine("Matrix dependency test: " + DepTest); Assert.LessOrEqual(DepTest, 0.0); // diagnostic output Console.WriteLine("Number of Agglomerations (all species): " + Agg.TotalNumberOfAgglomerations); Console.WriteLine("Number of Agglomerations (species 'B'): " + Agg.GetAgglomerator(LsTrk.GetSpeciesId("B")).AggInfo.SourceCells.NoOfItemsLocally.MPISum()); // operator auswerten: double[] x = new double[Affine.Length]; BLAS.daxpy(x.Length, 1.0, Affine, 1, x, 1); OperatorMatrix.SpMVpara(1.0, u.CoordinateVector, 1.0, x); MassInv.SpMV(1.0, x, 0.0, du_dx.CoordinateVector); Agg.GetAgglomerator(LsTrk.GetSpeciesId("B")).Extrapolate(du_dx.Mapping); // markieren, wo ueberhaupt A und B sind Bmarker.AccConstant(1.0, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); Amarker.AccConstant(+1.0, LsTrk.Regions.GetSpeciesSubGrid("A").VolumeMask); if (usePhi0 && usePhi1) { Xmarker.AccConstant(+1.0, LsTrk.Regions.GetSpeciesSubGrid("X").VolumeMask); } // compute error ERR.Clear(); ERR.Acc(1.0, du_dx_Exact, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); ERR.Acc(-1.0, du_dx, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); double L2Err = ERR.L2Norm(LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); Console.WriteLine("L2 Error: " + L2Err); XERR.Clear(); XERR.GetSpeciesShadowField("B").Acc(1.0, ERR, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); double xL2Err = XERR.L2Norm(); Console.WriteLine("L2 Error (in XDG space): " + xL2Err); // check error double ErrorThreshold = 1.0e-1; if (this.MomentFittingVariant == XQuadFactoryHelper.MomentFittingVariants.OneStepGaussAndStokes) { ErrorThreshold = 1.0e-6; // HMF is designed for such integrands and should perform close to machine accuracy; on general integrands, the precision is different. } bool IsPassed = ((L2Err <= ErrorThreshold || this.THRESHOLD <= ErrorThreshold) && xL2Err <= ErrorThreshold); if (IsPassed) { Console.WriteLine("Test PASSED"); } else { Console.WriteLine("Test FAILED: check errors."); //PlotCurrentState(phystime, TimestepNo, 3); } if (TimestepNo > 1) { if (this.THRESHOLD > ErrorThreshold) { // without agglomeration, the error in very tiny cut-cells may be large over the whole cell // However, the error in the XDG-space should be small under all circumstances Assert.LessOrEqual(L2Err, ErrorThreshold, "DG L2 error of computing du_dx"); } Assert.LessOrEqual(xL2Err, ErrorThreshold, "XDG L2 error of computing du_dx"); } // return/Ende base.NoOfTimesteps = 17; //base.NoOfTimesteps = 2; dt = 0.3; return(dt); }