/// <summary> /// Returns the source which contains Au_{n-1}+b for the Operator and the field provided in the constructor /// </summary> /// <param name="k">results of Ay+b</param> /// <param name="dt">optional scaling by time step size</param> public void ExplicitEulerSource(SubgridCoordinateMapping u, out double[] convectiveSource, double dt) { convectiveSource = new double[u.subgridCoordinates.Length]; u.Compress(); SubgridOperatorMatr.SpMVpara <double[], double[]>(dt, u.subgridCoordinates, 1.0, convectiveSource); BLAS.daxpy(SubgridAffine.Length, dt, SubgridAffine, 1, convectiveSource, 1); }
/// <summary> /// method for setting up the timestepper, i.e. the necessary /// </summary> protected void Setup1(Context ctx, ISparseSolver solver, bool[] temporalOp, MsrMatrix spatialOpMtx, IList <double> spatialOpAffine, SubGrid subgrid, SubgridCoordinateMapping fields, double InitialDeltat) { // check operator and arguments if (spatialOpMtx.NoOfRows != spatialOpMtx.NoOfCols) { throw new ArgumentException("matrix must be quadratic.", "spatialOpMtx"); } if (spatialOpMtx.NoOfRows != fields.GlobalCount) { throw new ArgumentException("matrix size must be equal to the GlobalCount of fields mapping", "fields,spatialOpMtx"); } if (spatialOpMtx.RowPartition.LocalLength != fields.NUpdate) { throw new ArgumentException("number of locally stored matrix rows nust be equal to NUpdate of fields mapping.", "fields,spatialOpMtx"); } if (spatialOpAffine.Count < fields.NUpdate) { throw new ArgumentException("length affine offset vector must be equal or larger than NUpdate of the mapping", "spatialOpAffine"); } m_Context = ctx; m_Solver = solver; m_Subgrid = subgrid; m_SubgridMapping = fields; m_AffineOffset1 = spatialOpAffine.ToArray(); this.temporalOp = temporalOp; BLAS.dscal(m_AffineOffset1.Length, -1.0, m_AffineOffset1, 1); { // check temporal operator // ----------------------- if (m_SubgridMapping.Fields.Count != temporalOp.Length) { throw new ArgumentException( "lenght of temporalOp must be equal to number of domain/codomain variables of the spatial differential operator", "temporalOp"); } m_SubgridMapping.Compress(); m_SubgridDGCoordinates = m_SubgridMapping.subgridCoordinates; m_SubgridMapping.SetupSubmatrix(m_AffineOffset1, spatialOpMtx, out m_CompressedAffine, out m_CompressedMatrix); bool timedep = false; bool fullyTimeDep = true; foreach (bool b in temporalOp) { timedep = timedep || b; fullyTimeDep = fullyTimeDep & b; } if (!timedep) { throw new ArgumentException("At least one equation must be time-dependent; one entry in temporalOp must be true;", "temporalOp"); } DefineMatrix(m_CompressedMatrix, InitialDeltat); } }
/// <summary> /// Constructor for an explicit Euler scheme operating on subgrids supporting parameters. /// </summary> /// <param name="ctx"></param> /// <param name="subgridMapping">Coordinate Mapping on the subgrid</param> /// <param name="Parameters">Optional parameters which have to match the matrix dimensions</param> /// <param name="operatorMatrix">Matrix of the differential operator</param> /// <param name="affine">Affine part of the operator matrix</param> public ExplicitEulerSubgrid(Context ctx, SubgridCoordinateMapping subgridMapping, CoordinateMapping Parameters, MsrMatrix operatorMatrix, double[] affine) { using (new ilPSP.Tracing.FuncTrace()) { m_Context = ctx; m_SubgridMapping = subgridMapping; m_DGCoordinates = subgridMapping.subgridCoordinates; m_Parameters = Parameters; IList <Field> ParameterFields = (m_Parameters == null) ? (new Field[0]) : m_Parameters.Fields; m_SubgridMapping.SetupSubmatrix(affine, operatorMatrix, out subgridAffine, out subgridMatrix); } }
public OperatorSplitting(Context context, SubGrid subgrid , SubgridCoordinateMapping v) { m_Context = context; sourcevariable = v; sourcevariable.Compress(); Partition part = new Partition(v.NUpdate, 1); affine = new double[v.NUpdate]; //opmatr= new MsrMatrix(part,v.MaxTotalNoOfCoordinatesPerCell * (int)m_Context.GridDat.GlobalNoOfCells); opmatr = new MsrMatrix(part); m_Subgrid = subgrid; }
public ExplicitConvection(Context context, string variablename, SubGrid subgrid, int basisdegree , SubgridCoordinateMapping v) { m_Context = context; name = variablename; convectionop = new SpatialOperator(1, 3, 1, name, "u", "v", "w"); convectionop.EquationComponents["Density"].Add(new SurfaceConvectionUpwinding(new string[] { "u", "v", "w" }, subgrid, new string[] { name }, 0)); convectionop.Commit(); CreepingFlowFactory fact = new CreepingFlowFactory(m_Context, basisdegree); fact.CreateFlowFields(out creepingFlow); Partition part = new Partition(v.NUpdate); double[] affine = new double[v.NUpdate]; MsrMatrix opmatr = new MsrMatrix(part, v.MaxTotalNoOfCoordinatesPerCell * (int)m_Context.GridDat.GlobalNoOfCells); convectionop.ComputeMatrixEx(m_Context, v, new Field[] { creepingFlow[0], creepingFlow[1], creepingFlow[2] }, v, opmatr, affine, false, subgrid); v.SetupSubmatrix(affine, opmatr, out SubgridAffine, out SubgridOperatorMatr); }
/// <summary> /// Constructor /// </summary> /// <param name="ctx">Context</param> /// <param name="subgridMapping">Coordinate Mapping that may be restricted to the subgrid</param> /// <param name="operatorMatrix">Matrix associated with the differential operator</param> /// <param name="affine">Affine part of the operator matrix</param> /// <param name="scheme">Runge-Kutta scheme, available in <see cref="RungeKutta"/></param> /// <param name="Parameters"> /// optional parameter fields, can be null if the spatial operator contains no parameters; /// must match the parameter field list of the spatial operator, see <see cref="BoSSS.Foundation.SpatialOperator.ParameterVar"/> /// </param> public RungeKuttaSubgrid(Context ctx, SubgridCoordinateMapping subgridMapping, CoordinateMapping Parameters, MsrMatrix operatorMatrix, double[] affine, RungeKutta.RungeKuttaScheme scheme) : base(ctx, subgridMapping, Parameters, operatorMatrix, affine) { m_RKscheme = scheme; }
/// <summary> /// sets the <see cref="SubgridCoordinateMapping"/> to <paramref name="m"/> /// The new mapping is compressed again to the narrow band only. /// </summary> /// <param name="m"></param> protected void SetMapping(SubgridCoordinateMapping m) { m_SubgridMapping = m; m_SubgridMapping.Compress(); m_SubgridDGCoordinates = m_SubgridMapping.subgridCoordinates; }
/// <summary> /// Factory for implicit time stepping schemes that operate on a subgrid /// </summary> /// <typeparam name="TimeStepperType"></typeparam> /// <param name="ctx"></param> /// <param name="solver">Solver</param> /// <param name="temporalOp">true for each component which should be time dependent</param> /// <param name="spatialOpMtx">Matrix associated with the differential operator</param> /// <param name="spatialOpAffine">Affine part of the differential operator</param> /// <param name="subgrid">Subgrid which the problem should be computed on</param> /// <param name="fields">Mapping on the subgrid</param> /// <param name="Initialdt">Time step size</param> /// <returns></returns> public static TimeStepperType Factory <TimeStepperType>(Context ctx, ISparseSolver solver, bool[] temporalOp, MsrMatrix spatialOpMtx, IList <double> spatialOpAffine, SubGrid subgrid, SubgridCoordinateMapping fields, double Initialdt) where TimeStepperType : ImplicitTimeStepperSubgrid, new() { TimeStepperType timestepper = new TimeStepperType(); timestepper.Setup1(ctx, solver, temporalOp, spatialOpMtx, spatialOpAffine, subgrid, fields, Initialdt); return(timestepper); }
/// <summary> /// Constructor for implicit time stepping scheme using a predefined discrization of the differential operator /// </summary> /// <param name="ctx"></param> /// <param name="solver">Solver</param> /// <param name="temporalOp">true for each component which should be time dependent</param> /// <param name="spatialOpMtx">Matrix associated with the differential operator</param> /// <param name="spatialOpAffine">Affine part of the differential operator</param> /// <param name="subgrid">Subgrid which the problem should be computed on</param> /// <param name="fields">Mapping on the subgrid</param> /// <param name="Initialdt">Time step size</param> public ImplicitTimeStepperSubgrid(Context ctx, ISparseSolver solver, bool[] temporalOp, MsrMatrix spatialOpMtx, IList <double> spatialOpAffine, SubGrid subgrid, SubgridCoordinateMapping fields, double Initialdt) { Setup1(ctx, solver, temporalOp, spatialOpMtx, spatialOpAffine, subgrid, fields, Initialdt); }
/// <summary> /// Constructor for an implicit Euler scheme that is operating on a subgrid. /// No ISparseExt solver is needed! /// </summary> /// <param name="ctx"></param> /// <param name="solver">The ISparse solver to be used</param> /// <param name="spatialOpMtx"> Matrix of the operator on the full domain</param> /// <param name="spatialOpAffine">Affine part of the operator defined on the full domain</param> /// <param name="subgrid">Subgrid where the computation is assumed to be performed</param> /// <param name="fields">Subgrid Mapping of all fields</param> public ImplicitEulerSubgrid(Context ctx, ISparseSolver solver, MsrMatrix spatialOpMtx, IList <double> spatialOpAffine, SubGrid subgrid, SubgridCoordinateMapping fields, double Initialdt) : this(ctx, solver, AllTrue(fields.Fields.Count), spatialOpMtx, spatialOpAffine, subgrid, fields, Initialdt) { }
/// <summary> /// Constructor for an explicit Euler scheme operating on subgrids(without parameters). /// </summary> /// <param name="ctx"></param> /// <param name="subgridMapping">Coordinate Mapping on the subgrid</param> /// <param name="operatorMatrix">Matrix of the differential operator</param> /// <param name="affine">Affine part of the operator matrix</param> public ExplicitEulerSubgrid(Context ctx, SubgridCoordinateMapping subgridMapping, MsrMatrix operatorMatrix, double[] affine) : this(ctx, subgridMapping, null, operatorMatrix, affine) { }
/// <summary> /// Evaluation of the operator on the subgrid by matrix-vector product. There might be a more efficient method.... /// </summary> /// <param name="k">results of Ay+b</param> /// <param name="dt">optional scaling by time step size</param> protected void Evaluate(SubgridCoordinateMapping u, double[] k, double dt) { SubgridOperatorMatr.SpMVpara <double[], double[]>(dt, u.subgridCoordinates, 1.0, k); BLAS.daxpy(SubgridAffine.Length, dt, SubgridAffine, 1, k, 1); }
public void ExplicitEulerSource(SubgridCoordinateMapping u, out double[] convectiveSource, double dt) { convectiveSource = new double[u.subgridCoordinates.Length]; Evaluate(u, convectiveSource, dt); }