Skip to content
Snippets Groups Projects
Commit 392a6c03 authored by Derek Steinmoeller's avatar Derek Steinmoeller
Browse files

Merge branch 'ds/NHSW' into 'master'

Ds/nhsw

See merge request !1
parents 9d7fa5dc b45bd1d2
No related branches found
No related tags found
1 merge request!1Ds/nhsw
......@@ -73,12 +73,21 @@ BASECASE_OBJS := $(addprefix BaseCase/,$(notdir $(BASECASE_CPPS:.cpp=.o)))
SCIENCE_CPPS := $(wildcard Science/*.cpp)
SCIENCE_OBJS := $(addprefix Science/,$(notdir $(SCIENCE_CPPS:.cpp=.o)))
# Non-hydrostatic shallow water deps
NH_BASE := TArray.o Parformer.o T_util.o gmres.o gmres_2d_solver.o gmres_1d_solver.o Splits.o Par_util.o gmres.o grad.o multigrid.o chebmat.o
##
## Recipes
##
NSIntegrator.o: NSIntegrator.cpp NSIntegrator_impl.cc
NHchannel_x: NHchannel.o $(NH_BASE)
$(LD) $(LDFLAGS) -o $@ $^ $(LDLIBS)
NHdonut2d_x: NHdonut2d.o $(NH_BASE)
$(LD) $(LDFLAGS) -o $@ $^ $(LDLIBS)
tests/test%.o: tests/test%.cpp
$(MPICXX) $(CFLAGS) -o $@ -c $<
......
This diff is collapsed.
#include <blitz/array.h>
#include <mpi.h>
#include "T_util.hpp"
#include "Par_util.hpp"
#include "TArray.hpp"
#include <stdio.h>
#include <math.h>
#include <iostream>
#include <vector>
#include "Splits.hpp"
#include <random/normal.h>
#include <fstream>
#include "Parformer.hpp"
#include <stdlib.h>
#include <stdio.h>
#include "gmres.hpp"
#include "gmres_2d_solver.hpp"
#include "gmres_2d_solver_impl.hpp"
#include "grad.hpp"
#include "chebmat.hpp"
#include "lserk4.hpp"
using namespace TArrayn;
using namespace Transformer;
using blitz::Array;
using blitz::TinyVector;
using blitz::GeneralArrayStorage;
using blitz::Range;
using ranlib::Normal;
using namespace std;
// Grid size
#define N_R 32
#define N_THETA 128
// Grid lengths
#define R_OUTER 8000.0
#define R_INNER 1000.0
#define L_R (R_OUTER - R_INNER)
// Defines for physical parameters
#define G (0.03*9.81)
#define EARTH_OMEGA (2*M_PI/(24*3600))
#define EARTH_RADIUS (6371e3)
#define LATITUDE (M_PI/2)
#define ROT_F (7.8828e-5)
//#define ROT_F (2*EARTH_OMEGA*sin(LATITUDE))
#define ROT_B (0*2*EARTH_OMEGA*cos(LATITUDE)/EARTH_RADIUS)
#define H0 (12.8)
#define H_DEPTH (H0) //const for now.
// Timestep parameters
#define FINAL_T 1200.0
//#define FINAL_T (3600.0*24.0*3.0)
#define INITIAL_T 0.0
#define SAFETY_FACTOR 0.25
#define NUMOUTS 200.0
// Filtering parameters
#define FILTER_ON true
#define F_CUTOFF 0.4 //0.4 usually
#define F_ORDER 4
#define F_STREN (-0.33)
// GMRES parameters
#define MAXIT 30
#define TOL 1.0e-8
#define RESTARTS 1
#define NOISYGMRES true
// do we want to output at all?
#define OUTFLAG true
#define NHOUTFLAG true //output auxiliary variable, z?
// Blitz index placeholders
blitz::firstIndex ii;
blitz::secondIndex jj;
blitz::thirdIndex kk;
int main(int argc, char ** argv) {
// Initialize MPI
MPI_Init(&argc, &argv);
// make grids
Array<double,1> theta(N_THETA), r(N_R);
// matrices for dealing with imposing Neuman bc's on eta
Array<double,2> Dr(N_R,N_R);
Array<double,2> neumatinv(2,2);
double dtheta = 2*M_PI/N_THETA;
theta = (ii+0.5)*dtheta;
r = -cos(M_PI*ii/(N_R-1));
// get cheb matrix on [-1,1]
Dr = chebmat(N_R);
//adjust from standard grids to physical grids
r = L_R*((r(ii)+1.0)/2.0) + R_INNER;
//apply Jacobian to get cheb matrix on physical grid.
Dr = (-2.0/L_R)*Dr;
// put entries into 2x2 Neuman-imposing matrix
double neudet = Dr(0,0)*Dr(N_R-1,N_R-1)-Dr(0,N_R-1)*Dr(N_R-1,0);
neumatinv(0,0) = Dr(N_R-1,N_R-1);
neumatinv(0,1) = -Dr(0,N_R-1);
neumatinv(1,0) = -Dr(N_R-1,0);
neumatinv(1,1) = Dr(0,0);
neumatinv = (1/neudet)*neumatinv;
// Get parameters for local array storage
TinyVector<int,3> local_lbound, local_extent;
GeneralArrayStorage<3> local_storage;
local_lbound = alloc_lbound(N_THETA,1,N_R);
local_extent = alloc_extent(N_THETA,1,N_R);
local_storage = alloc_storage(N_THETA,1,N_R);
// Necessary FFT transformers (for filtering)
TransWrapper XY_xform(N_THETA,1,N_R,FOURIER,NONE,CHEBY);
// Allocate space for flow variables
DTArray utheta(local_lbound,local_extent,local_storage);
DTArray ur(local_lbound,local_extent,local_storage);
DTArray eta(local_lbound,local_extent,local_storage);
// Runge-Kutta residual storage.
DTArray res_utheta(local_lbound,local_extent,local_storage);
DTArray res_ur(local_lbound,local_extent,local_storage);
DTArray res_eta(local_lbound,local_extent,local_storage);
// Allocate array for depth profile
DTArray H(local_lbound,local_extent,local_storage);
// Common geometric factor (1/r)
DTArray rinv(local_lbound,local_extent,local_storage), r2d(local_lbound,local_extent,local_storage);
rinv = 1/r(kk);
r2d = r(kk);
// set equal to expression defined in pre-processing
H = H_DEPTH;
double Hmax = max(H);
// set initial conditions
//etaq = (H0/4)*exp(-.1*((y(kk)/1e2)*(y(kk)/1e2))
// -1*((x(ii)-0.5*Lx)/3e2)*((x(ii)-0.5*Lx)/3e2));
eta = (H0/4.0)*(r2d/R_OUTER)*(cos(theta(ii)) + 0*kk);
ur = 0.0;
utheta = 0.0;
int tstep = 0; //time-step counter
// Compute time-stepping details
double etamax = psmax(max(abs(eta)));
double c0max = sqrt(G*(Hmax + etamax));
double dr = fabs(r(1)-r(0));
double dt = SAFETY_FACTOR*fmin(dr,R_INNER*dtheta)/c0max;
// this assumes uniform spacing in Cheby (y-)direction, should fix.
double t = INITIAL_T;
double numsteps = (FINAL_T - INITIAL_T)/dt;
int outputinterval = (int) floor(numsteps/NUMOUTS);
// don't allow zero.
outputinterval = max(outputinterval,1);
if (master()) printf("Using timestep of %gs, with final time of %gs\n",dt,FINAL_T);
if (master()) printf("Reference: C_0_max = %g m/s, dr = %g m, d(arc) = %g m\n",c0max,L_R/N_R,R_INNER*dtheta);
// initialize output stream for times file, and write out initial time.
ofstream timefile;
if (master()) {
timefile.open("times", ios::out | ios::trunc);
timefile << t << "\n";
timefile.close();
}
// Allocate space for auxiliary elliptic variable
DTArray z(local_lbound,local_extent,local_storage);
if (OUTFLAG == true) {
write_reader(eta, "eta", true);
write_reader(ur, "ur", true);
write_reader(utheta, "utheta", true);
write_reader(H, "H");
write_array(eta,"eta",tstep/outputinterval);
write_array(ur,"ur",tstep/outputinterval);
write_array(utheta,"utheta",tstep/outputinterval);
write_array(H,"H");
if (NHOUTFLAG == true)
write_reader(z,"z",true);
}
// Allocate space for some temporaries
DTArray temp1(local_lbound,local_extent,local_storage),
temp2(local_lbound,local_extent,local_storage),
temp3(local_lbound,local_extent,local_storage),
temp4(local_lbound,local_extent,local_storage),
temp5(local_lbound,local_extent,local_storage);
// Allocate space for forcing vector
DTArray Fx(local_lbound,local_extent,local_storage),
Fy(local_lbound,local_extent,local_storage);
// Allocate space for elliptic problem variable coefficients.
DTArray cx(local_lbound,local_extent,local_storage),
cxx(local_lbound,local_extent,local_storage),
cy(local_lbound,local_extent,local_storage),
cyy(local_lbound,local_extent,local_storage);
// write grids to file
temp1 = theta(ii) + 0.0*kk;
write_array(temp1,"thetagrid");
write_reader(temp1,"thetagrid");
write_array(r2d,"rgrid");
write_reader(r2d,"rgrid");
// Build gradient operator, calling seq: szx,szy,szz, s_exp in x,y,z
Grad mygrad(N_THETA,1,N_R,FOURIER,NONE,CHEBY);
// set (constant) jacobian based on our rectangular grid.
mygrad.set_jac(firstDim,firstDim,2*M_PI);
mygrad.set_jac(secondDim,secondDim,1.0);
mygrad.set_jac(thirdDim,thirdDim,-2/L_R); //Cheby grids need the '-'.
// set coefficients of elliptic problem
cyy = (H*H)/6; //coeff of z_rr
cxx = (rinv*rinv)*cyy; //coeff of z_{\theta\theta}
//cx & cy (coefs of z_\theta and z_r) are calculated numerically.
mygrad.setup_array(&cxx, FOURIER,NONE,CHEBY);
mygrad.get_dx(&cx); //get theta derivative, store in cx. This is d(H2o6/r^2)/d\theta
temp1 = r2d*cyy;
mygrad.setup_array(&temp1, FOURIER,NONE,CHEBY);
mygrad.get_dz(&temp2); //get r derivative, store in temp2. This is d(r*H2o6)/dr
cy = rinv*temp2; // needs to be divided by r to complete the Christoffel symbol
// create gmres-2d solver interface
Cheb_2dmg a(N_THETA, N_R);
int itercount; //gmres iteration count
double dbc=0.0;
double nbc=1.0;
// set variable coefs (must do before set_grad)
a.set_ci(&cx,1);
a.set_ci(&cxx,2);
a.set_ci(&cy,3);
a.set_ci(&cyy,4);
a.set_grad(&mygrad); // need to do this before set_bc
// set_bc, calling sequence is helm, dir., neu., S_EXP ('x' series expans)
a.set_bc(-1.0,dbc,nbc,FOURIER,dbc,nbc);
// Initialize structures that are used by gmres object
fbox * init_r = a.alloc_resid();
ubox * final_u = a.alloc_basis();
//Build the solver object with template class
GMRES_Solver<Cheb_2dmg> mysolver(&a);
// set forcing
Fx = 0;
Fy = 0;
double starttime = MPI_Wtime();
// Runge-Kutta steps
if (master()) printf("Entering main loop... \n");
while (t < FINAL_T) {
//make nice references for arrays used
DTArray &temp = temp1, &rhs=temp2, &a1 = temp3, &a2 = temp4, &dump = temp5;
//update forcing:
Fx = 0.0; Fy = 0.0;
// Initialize with zero RK residual.
res_utheta = 0.0; res_ur = 0.0; res_eta = 0.0;
for (int i=0; i < LSERK4::numStages; i++) {
//step eta:
//eta_p = eta_m - dt*((eta_n u)_x + (eta_n v)_y)
temp = -r2d*((H+eta)*ur);
mygrad.setup_array(&temp,FOURIER,NONE,CHEBY);
mygrad.get_dz(&rhs, false); //get r derivative, store in rhs
rhs = rinv*rhs; // (1/r)*d(r*F_r)/dr where F_r = (H+eta)u_r
temp = -rinv*(H+eta)*utheta; // d((1/r)*(H+eta)*utheta)/dtheta
mygrad.setup_array(&temp,FOURIER,NONE,CHEBY);
mygrad.get_dx(&rhs, true); //get theta derivative. This time, we sum with old result (cumulative = true).
res_eta = LSERK4::rk4a[i]*res_eta + dt*rhs;
eta += LSERK4::rk4b[i]*res_eta;
// filter free surface
if (FILTER_ON == true)
filter3(eta,XY_xform,FOURIER,NONE,CHEBY,F_CUTOFF,F_ORDER,F_STREN);
//impose explicit Neuman conditions on eta.
double neurhs1 = 0.0;
double neurhs2 = 0.0;
for (int j=local_lbound(1); j<(local_lbound(1)+local_extent(1)); j++) {
neurhs1 = sum(-Dr(0,Range(1,N_R-2))*eta(j,0,Range(1,N_R-2)));
neurhs2 = sum(-Dr(N_R-1,Range(1,N_R-2))*eta(j,0,Range(1,N_R-2)));
// perform 2x2 matrix product at each point
eta(j,0,0) = neumatinv(0,0)*neurhs1 + neumatinv(0,1)*neurhs2;
eta(j,0,N_R-1) = neumatinv(1,0)*neurhs1+ neumatinv(1,1)*neurhs2;
}
//form a1 and a2:
//nonlinear advection
mygrad.setup_array(&ur,FOURIER,NONE,CHEBY);
mygrad.get_dz(&temp,false); //d(ur)/dr
a1 = -ur*temp; // ur * d(ur)/dr
mygrad.get_dx(&temp,false); //d(ur)/d\theta
a1 -= rinv*utheta*temp; // (1/r)*utheta * d(ur)/dtheta
a1 += rinv*utheta*utheta; // Geometric term.
// a2
mygrad.setup_array(&utheta,FOURIER,NONE,CHEBY);
mygrad.get_dz(&temp,false); //d(utheta)/dr
a2 = -ur*temp; // ur * d(utheta)/dr
mygrad.get_dx(&temp,false); //d(utheta)/d\theta
a2 -= rinv*utheta*temp; // (1/r)(utheta) * d(utheta)/dtheta
a2 -= rinv*utheta*ur; // Geometric term.
// pressure gradient
mygrad.setup_array(&eta,FOURIER,NONE,CHEBY);
mygrad.get_dz(&temp,false); // d(eta)/dr
a1 -= G*temp;
mygrad.get_dx(&temp,false); // (1/r)*d(eta)_/dtheta
a2 -= G*rinv*temp;
// Coriolis and forcing
a1 += (ROT_F+ROT_B*r(kk)*sin(theta(ii)))*utheta + Fx;
a2 -= (ROT_F+ROT_B*r(kk)*sin(theta(ii)))*ur + Fy;
//done forming a1 & a2
// Perform Runge-Kutta incremental step..
res_ur = LSERK4::rk4a[i]*res_ur + dt*a1;
res_utheta = LSERK4::rk4a[i]*res_utheta + dt*a2;
ur += LSERK4::rk4b[i]*res_ur;
utheta += LSERK4::rk4b[i]*res_utheta;
} // RK-stages loop.
//compute - divergence of \vec{a}, store in rhs for helmholtz problem
dump = r2d*a1;
mygrad.setup_array(&dump,FOURIER,NONE,CHEBY);
mygrad.get_dz(&temp,false);
temp = rinv*temp; // (1/r)*d(r*a1)/dr
rhs = -temp;
mygrad.setup_array(&a2,FOURIER,NONE,CHEBY);
mygrad.get_dx(&temp,false);
rhs -= rinv*temp; // (1/r)*d(a2)/dtheta
//set bc's on RHS
rhs(Range::all(),0,0) = -a1(Range::all(),0,0) / cyy(Range::all(),0,0);
rhs(Range::all(),0,N_R-1) = -a1(Range::all(),0,N_R-1) / cyy(Range::all(),0,N_R-1);
// check for blow-up
double sanity = pvmax(rhs);
if (master()) {
if (isnan(sanity) || isinf(sanity) || sanity > 1.e10) {
cout << "sanity: " << sanity << endl;
cout << "Numerical blow-up has occurred, gg." << endl;
return 1;
}
}
//solve for z (non-hydrostatic pressure) with gmres
*init_r->gridbox = rhs;
itercount = mysolver.Solve(final_u,init_r,TOL,MAXIT,RESTARTS);
if (master() && itercount < 0) {
cout << "GMRES did not converged, with status: " << itercount << endl;
cout << "Terminating..." << endl;
return 1;
}
if (master() && NOISYGMRES == true)
cout << "GMRES converged after " << itercount << " iterations." << endl;
z = *final_u->gridbox;
//compute gradient of NH pressure, then time-step it & velocities
mygrad.setup_array(&z,FOURIER,NONE,CHEBY);
mygrad.get_dz(&temp,false); // dz/dr
ur += dt*(a1 + cyy*temp); //cyy = H2o6
mygrad.get_dx(&temp,false); // (1/r)dz/d\theta
utheta += dt*(a2 + cyy*rinv*temp); //cyy = H2o6
//filter velocity
if (FILTER_ON == true) {
filter3(ur, XY_xform,FOURIER,NONE,CHEBY,F_CUTOFF,F_ORDER,F_STREN);
filter3(utheta, XY_xform,FOURIER,NONE,CHEBY,F_CUTOFF,F_ORDER,F_STREN);
}
//Impose BCs on v at r=R_INNER,R_OUTER (Dirichlet -> no normal flow)
ur(Range::all(), 0, 0) = 0;
ur(Range::all(), 0, N_R-1) = 0;
//increment time.
tstep++;
t+=dt;
// Check if it's time to output
if(!(tstep % outputinterval) && OUTFLAG == true){
if (master())
cout << "outputting at t=" << t << "\n";
write_array(ur ,"ur", tstep/outputinterval);
write_array(utheta,"utheta", tstep/outputinterval);
write_array(eta ,"eta", tstep/outputinterval);
if (NHOUTFLAG == true)
write_array(z,"z",tstep/outputinterval);
//append current time to times file
if (master()) {
timefile.open ("times", ios::out | ios::app);
timefile << t << "\n";
timefile.close();
}
}
//generic log-style text output to screen
if(!(tstep % 100) || tstep < 20) {
if (master())
printf("Completed time %g (tstep %d)\n",t,tstep);
double mur = psmax(max(abs(ur))), mutheta = psmax(max(abs(utheta))), meta = psmax(max(abs(eta)));
if (master())
printf("Max ur %g, utheta %g, eta %g\n",mur,mutheta,meta);
}
} //end time-stepping loop
double endtime = MPI_Wtime();
double runtime = endtime-starttime;
if (master()) {
printf("Finished run. Run-time = %g s \n", runtime);
timefile.open ("runtime", ios::out | ios::trunc);
timefile << runtime << "\n";
timefile.close();
}
// free memory of residual and basis of gmres solver
a.free_resid(init_r);
a.free_basis(final_u);
MPI_Finalize();
return 0;
}
\ No newline at end of file
#include "chebmat.hpp"
using blitz::Array;
using blitz::TinyVector;
using blitz::GeneralArrayStorage;
using blitz::Range;
using blitz::firstDim;
using blitz::secondDim;
// Blitz index placeholders
//blitz::firstIndex ii;
//blitz::secondIndex jj;
//blitz::thirdIndex kk;
Array<double,2> chebmat(int Nx) {
blitz::firstIndex ii;
blitz::secondIndex jj;
Array<double,1> x(Nx);
x = cos(M_PI*ii/(Nx-1));
Array<double,2> myI(Nx,Nx);
Array<double,2> Dx(Nx,Nx);
Array<double,2> myX(Nx,Nx);
Array<double,1> myc(Nx);
Array<double,2> mydX(Nx,Nx);
Array<double,2> mysum(Nx,Nx);
// build identity matrix (any way to clean this up?)
for (int j=0; j<Nx; j++)
myI(j,j) = 1.0;
myc(0) = 2.0;
myc(Range(1,Nx-2)) = 1.0;
myc(Nx-1) = 2.0;
myc = myc*pow(-1.0,ii);
// like the "repmat" step - clean this up with blitz functionality?
for (int j=0; j<Nx; j++)
myX(Range::all(),j) = x;
mydX = myX-myX.transpose(secondDim,firstDim);
Dx = myc(ii)*(1/myc(jj)); //outer product
Dx = Dx/(mydX + myI);
// build diagonal matrix of row-sums
// perhaps try: sum(Dx(i,Range::all()))
double rowsum;
for (int i=0; i<Nx; i++) {
rowsum = 0.0;
for (int j=0; j<Nx; j++) {
rowsum+= Dx(i,j);
}
mysum(i,i) = rowsum;
}
Dx = Dx - mysum;
return Dx;
}
\ No newline at end of file
#pragma once
#include <blitz/array.h>
#include <math.h>
blitz::Array <double,2> chebmat(int Nx);
\ No newline at end of file
......@@ -249,6 +249,11 @@ void Cheb_2dmg::build_precond() {
// Now, get the Jacobian from grad_3d, carve off a relevant 2D slice, and
// assign it to the 2d gradient's jacobian.
double c_jac; // entry for constant jacobian
double myJ11=0.0, myJ22=0.0; //Derek: placeholders for constant jacobians.
// We assume these two are consts, and
// J12=J21= 0
// xx
grad_3d->get_jac(firstDim,firstDim,c_jac,v_jac);
if (!v_jac) {// constant
......@@ -358,6 +363,7 @@ void Cheb_2dmg::build_precond() {
/* J11 */
grad_2d->get_jac(firstDim,firstDim,c_jac,v_jac);
if (!v_jac) { // constant
myJ11 = c_jac; //Derek
t1 = c_jac;
t4 = 0;
} else {
......@@ -368,6 +374,7 @@ void Cheb_2dmg::build_precond() {
// J22
grad_2d->get_jac(thirdDim,thirdDim,c_jac,v_jac);
if (!v_jac) { // constant
myJ22 = c_jac; //Derek
t2 = c_jac;
t5 = 0;
} else {
......@@ -403,6 +410,18 @@ void Cheb_2dmg::build_precond() {
t2 = t2*t2 + (*v_jac)*(*v_jac);
}
//Derek: try generalizing it for any box-domain:
DTArray& helmCoef1 = *c1;
DTArray& helmCoef2 = *c2;
DTArray& helmCoef3 = *c3;
DTArray& helmCoef4 = *c4;
t1 = myJ11*myJ11*helmCoef2;
t2 = myJ22*myJ22*helmCoef4;
t4 = myJ11*helmCoef1;
t5 = myJ22*helmCoef3;
// Now, t1->t5 are set properly, so assign them to the MG operator
// fprintf(stderr,"PROBLEM SETUP\n");
// cerr << t1_2d << t2_2d << t3_2d << t4_2d << t5_2d;
......@@ -820,6 +839,23 @@ void Cheb_2dmg::set_bc(double h, double zdir, double zneum, S_EXP sym_type, doub
mg_precond->bc_setup(secondDim,dir_z,db_bot,da_bot,
dir_z,db_top,da_top);
}
// Derek method for setting the object's variable coefs. (all pointers)
void Cheb_2dmg::set_ci(DTArray * cin, int i) {
if (i==1)
c1 = cin;
else if (i==2)
c2 = cin;
else if (i==3)
c3 = cin;
else if (i==4)
c4 = cin;
else
fprintf(stderr,"set_ci\n");
return;
}
void Cheb_2dmg::precondition(fbox * & rhs, ubox * & soln) {
......@@ -930,6 +966,9 @@ void Cheb_2dmg::matrix_multiply(ubox * & lhs, fbox * & rhs) {
// return;
DTArray & dx = *temp1, & dz = *temp2;
//Derek: temporary for 2nd derivative.
DTArray& dd = *temp3;
Range all = Range::all();
// fprintf(stderr,"Matrix multiply\n");
......@@ -943,6 +982,17 @@ void Cheb_2dmg::matrix_multiply(ubox * & lhs, fbox * & rhs) {
grad_2d->get_dx(&dx);
grad_2d->get_dz(&dz);
DTArray& helmCoef1 = *c1;
DTArray& helmCoef2 = *c2;
DTArray& helmCoef3 = *c3;
DTArray& helmCoef4 = *c4;
// Derek: add first derivative terms to RHS
*rhs->gridbox = helmCoef1*dx + helmCoef3*dz;
//*rhs->gridbox = (*c1)*dx - (*c3)*dz; //Derek: Need minus if not
//doing setjac(-1,...) on
// z-coordinate.
// fprintf(stderr,"Dx\n");
// fprintf(stderr,"Dz\n");
......@@ -953,8 +1003,11 @@ void Cheb_2dmg::matrix_multiply(ubox * & lhs, fbox * & rhs) {
grad_2d->setup_array(&dx,SINE,NONE,CHEBY);
else
grad_2d->setup_array(&dx,type_x,NONE,CHEBY);
grad_2d->get_dx(rhs->gridbox);
// Derek:
grad_2d->get_dx(&dd);
*rhs->gridbox = *rhs->gridbox + helmCoef2*dd;
// lhs_zz
if (type_x == SINE)
grad_2d->setup_array(&dz,COSINE,NONE,CHEBY);
......@@ -962,7 +1015,10 @@ void Cheb_2dmg::matrix_multiply(ubox * & lhs, fbox * & rhs) {
grad_2d->setup_array(&dz,SINE,NONE,CHEBY);
else
grad_2d->setup_array(&dz,type_x,NONE,CHEBY);
grad_2d->get_dz(rhs->gridbox,true);
// Derek:
grad_2d->get_dz(&dd);
*rhs->gridbox = *rhs->gridbox + helmCoef4*dd;
// apply helmholtz
......
......@@ -81,6 +81,10 @@ class Cheb_2dmg : public GMRES_Interface < fbox * , ubox *> {
// used
DTArray * temp1, * temp2, * temp3;
//Derek: variable coefficients for more general problem
// c1 u_x + c2 u_xx + c3 u_z + c4 u_zz - u = f
DTArray * c1, * c2, * c3, * c4;
// C-ordered arrays for interfacing with the multigrid preconditioner
Array<double,2> * r2d, * s2d;
......@@ -137,6 +141,9 @@ class Cheb_2dmg : public GMRES_Interface < fbox * , ubox *> {
bool noisy() const;
// Derek: method for setting variable coefficients
void set_ci(DTArray * cin, int i);
private:
// Allocation data
vector<fbox *> r_free;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment