[petsc-users] Integrate PETSC with existing Fortran

Matthew Knepley knepley at gmail.com
Wed Jan 23 12:30:04 CST 2019


On Wed, Jan 23, 2019 at 1:27 PM Yaxiong Chen <chen2018 at purdue.edu> wrote:

> Hi Matt,
>
>
> I tried to modify the structure of present code and use KSP in the main
> program(optimal_mechanical_part).  Now the makefile is as following :
>
>
> #================================================================================
>
> OUT_DIR := HiDACX/
>
> ## Define mode
> ifndef MODE
>
> .PHONY: debug release
>
> debug: export MODE := debug
> release: export MODE := release
> debug release:
> @$(MAKE) $@
>
> clean:
> rm -rf $(OUT_DIR)obj/ $(OUT_DIR)lib/ $(OUT_DIR)bin/
>
> else
>
> ## SRCS = $(shell find ../../MyGit/Hidac_Fortran_Lib/ "*.f90")
> SRCS_DIR :=
> $(HOME)/Documents/HiDAC_Code_Git/HiDAC_Tao_NURBS_Petsc/Hidac_Fortran_Lib/
> FSRCS := \
> HiDAC_Utils/HiDAC_parameter.f90 \
> HiDAC_Utils/HiDAC_dataStructure.f90 \
> HiDAC_Utils/HiDAC_Math/HiDAC_mathFunction.f90 \
> HiDAC_Utils/HiDAC_Math/HiDAC_mathType.f90 \
> HiDAC_Utils/HiDAC_Math/HiDAC_tensor.f90 \
> HiDAC_Utils/HiDAC_Math/HiDAC_Math.f90 \
> HiDAC_Utils/HiDAC_Port/HiDAC_IGES.f90 \
> HiDAC_Utils/HiDAC_Port/IR_Precision.f90 \
> HiDAC_Utils/HiDAC_Port/Lib_Pack_Data.f90 \
> HiDAC_Utils/HiDAC_Port/Lib_Base64.f90 \
> HiDAC_Utils/HiDAC_Port/Lib_VTK_IO.f90 \
> HiDAC_Utils/HiDAC_Port/HiDAC_VTK.f90 \
> HiDAC_Utils/HiDAC_Port/HiDAC_Port.f90 \
> HiDAC_Utils/HiDAC_DistanceField/HiDAC_QhullInterface_F.f90 \
> HiDAC_Utils/HiDAC_DistanceField/HiDAC_BezierHull.f90 \
> HiDAC_Utils/HiDAC_DistanceField/HiDAC_DistanceField.f90 \
> HiDAC_Utils/HiDAC_Utils.f90 \
> HiDAC_Material/HiDAC_baseMat.f90 \
> HiDAC_Material/HiDAC_isoElasticMat.f90 \
> HiDAC_Material/HiDAC_thermoMat.f90 \
> HiDAC_Material/HiDAC_Material.f90 \
> HiDAC_Approximation/HiDAC_node.f90 \
> HiDAC_Approximation/HiDAC_nodeNet.f90 \
> HiDAC_Approximation/HiDAC_baseApproximation.f90 \
> HiDAC_Approximation/HiDAC_Approximation.f90 \
> HiDAC_Approximation/HiDAC_Approximation_NURBS/HiDAC_knot.f90 \
> HiDAC_Approximation/HiDAC_Approximation_NURBS/HiDAC_BezierApproximation.f90
> \
> HiDAC_Approximation/HiDAC_Approximation_NURBS/HiDAC_NURBSApproximation.f90
> \
> HiDAC_Approximation/HiDAC_Approximation_NURBS/HiDAC_THBApproximation.f90 \
> HiDAC_Domain/HiDAC_primitive.f90 \
> HiDAC_Domain/HiDAC_enrichment.f90 \
> HiDAC_Domain/HiDAC_composedDomain.f90 \
> HiDAC_Domain/HiDAC_Domain.f90 \
> HiDAC_Quadrature/HiDAC_quadType.f90 \
> HiDAC_Quadrature/HiDAC_directQuadrature.f90 \
> HiDAC_Quadrature/HiDAC_adaptiveQuadrature.f90 \
> HiDAC_Quadrature/HiDAC_Quadrature.f90 \
> HiDAC_Analysis/HiDAC_BC.f90 \
> HiDAC_Analysis/HiDAC_baseAnalysis.f90 \
> HiDAC_Analysis/HiDAC_elastostatics.f90 \
> HiDAC_Analysis/HiDAC_heatConduction.f90 \
> HiDAC_Analysis/HiDAC_phaseTransition.f90 \
> HiDAC_Analysis/HiDAC_fracture.f90 \
> HiDAC_Analysis/HiDAC_lineenrich.f90 \
> HiDAC_Analysis/HiDAC_sensitivityAnalysis.f90 \
> HiDAC_Analysis/HiDAC_Analysis.f90 \
> HiDAC_Solver/HiDAC_baseSolver.f90 \
> HiDAC_Solver/HiDAC_LapackSolver.f90 \
> HiDAC_Solver/HiDAC_MumpsSolver.f90 \
> HiDAC_Solver/HiDAC_Solver.f90 \
> IGES/read_iges_file.f90 \
> IGES/nrb.f90 \
> IGES/parse.f90 \
> IGES/map.f90 \
> DXF/dxf_reader.f90 \
>
> TSRCS_DIR := $(SRCS_DIR)HiDAC_Test/
>
> CSRCS := \
> HiDAC_Utils/HiDAC_DistanceField/HiDAC_QhullInterface_C.c \
>
> SRCS := $(addprefix $(SRCS_DIR), $(FSRCS) $(CSRCS))
>
> OBJS_DIR := $(OUT_DIR)obj/$(MODE)/
> LIB_DIR := $(OUT_DIR)lib/$(MODE)/
> EXE_DIR := $(OUT_DIR)bin/$(MODE)/
> EXT_DIR := /usr/local/
>
> #petsc lib
> include ${PETSC_DIR}/lib/petsc/conf/variables
> include ${PETSC_DIR}/lib/petsc/conf/rules
>
>
> #========================
>
> OBJS := $(addprefix $(OBJS_DIR), \
> $(patsubst %.f90, %.o, $(notdir $(FSRCS))) \
> $(patsubst %.c, %.o, $(notdir $(CSRCS))))
>
> LIB = $(LIB_DIR)libhidac.a
>
> VPATH := $(sort $(dir $(SRCS))) $(TSRCS_DIR)
>
> FC := mpif90
> CC := mpicc
> AR := ar rcs
>
> FCFLAGS := -cpp -J$(OBJS_DIR) # ${PETSC_KSP_LIB}# -fprofile-arcs
> -ftest-coverage
> CCFLAGS := -Dqh_QHpointer -ansi
> LFLAGS := -Wl,-rpath,$(EXT_DIR)lib, ${PETSC_KSP_LIB}
>
> FIDIR := \
> -I$(EXT_DIR)include \
> -I$(PETSC_DIR)/include
> CIDIR := \
> -I$(EXT_DIR)include/libqhull
>
> # sequence does matter
> EXTLIBS := \
> -L$(EXT_DIR)lib \
> -ldmumps -lmumps_common \
> -lparmetis -lmetis -lpord -lscalapack \
> -llapack -lblas -lqhullstatic_p \
> -lpetsc
>
>
> .PHONY: debug release all clean
>
> # Handle different mode using target-specific variable
> debug: FCFLAGS += -Wall -fcheck=bounds
> debug: CCFLAGS += -Wall
> debug: LFLAGS += -Wall -fcheck=bounds
>
> release: FCFLAGS += -O2
> release: CCFLAGS += -O2
> release: LFLAGS += -O2
>
> debug release: $(EXE_DIR)optimal_mechanical_part.hdc
>
> ## Rules
> .SUFFIXES: .hdc
> # Note the sequence of $(OBJS_DIR)%.o $(LIB) and $(EXTLIBS) does matter!
> Lib later!
> $(EXE_DIR)%.hdc: $(OBJS_DIR)%.o $(LIB)
> @mkdir -p $(EXE_DIR)
> echo $(OBJS_DIR)
> $(FC) $(LFLAGS) -o $@ $^ $(EXTLIBS) ${PETSC_KSP_LIB}
>
> $(LIB): $(OBJS)
> echo $(OBJS)
> @mkdir -p $(LIB_DIR)
> $(AR) $@ $^
>
> $(OBJS_DIR)%.o:
> @mkdir -p $(OBJS_DIR)
> @echo $@
> $(if $(findstring .f90, $<), \
> $(FC) $(FCFLAGS) $(FIDIR) -c $< -o $@, \
> $(CC) $(CCFLAGS) $(CIDIR) -c $< -o $@)
>
> ## Dependencies of files
> $(OBJS_DIR)optimal_mechanical_part.o: \
> optimal_mechanical_part.f90 \
> $(OBJS_DIR)read_iges_file.o \
> $(OBJS_DIR)HiDAC_Approximation.o \
> $(OBJS_DIR)HiDAC_Domain.o \
> $(OBJS_DIR)HiDAC_Material.o \
> $(OBJS_DIR)HiDAC_Quadrature.o \
> $(OBJS_DIR)HiDAC_Analysis.o \
> $(OBJS_DIR)HiDAC_Utils.o
>
> $(OBJS_DIR)read_iges_file.o: \
> read_iges_file.f90 \
> $(OBJS_DIR)HiDAC_Approximation.o \
> $(OBJS_DIR)HiDAC_Utils.o \
> $(OBJS_DIR)nrb.o \
> $(OBJS_DIR)parse.o \
> $(OBJS_DIR)map.o \
> $(OBJS_DIR)dxf_reader.o
>
> $(OBJS_DIR)test_main.o: \
> test_main.f90 \
> $(OBJS_DIR)read_iges_file.o
>
> $(OBJS_DIR)dxf_reader.o: \
> dxf_reader.f90 \
> $(OBJS_DIR)nrb.o
>
> $(OBJS_DIR)parse.o: \
> parse.f90 \
> $(OBJS_DIR)nrb.o
>
> $(OBJS_DIR)map.o: \
> map.f90 \
> $(OBJS_DIR)nrb.o
>
> <Other dependencies >
> # $(OBJS_DIR)HiDAC_PetscSolver.o : HiDAC_PetscSolver.f90
> $(OBJS_DIR)HiDAC_baseSolver.o $(OBJS_DIR)HiDAC_Utils.o
> # mpif90 -cpp -JHiDACX/obj/debug/ -Wall -ffree-line-length-0
> -Wno-unused-dummy-argument -g
> -I/Users/yaxiong/Downloads/petsc-3.9.4/include \
> # -I/Users/yaxiong/Downloads/petsc-3.9.4/arch-darwin-c-debug/include \
> # -I/opt/X11/include -o $(OBJS_DIR)HiDAC_PetscSolver.o \
> # /Users/yaxiong/Documents/HiDAC_Code_Git/HiDAC_Tao_NURBS_PETSC/Hidac_Fortran_Lib/HiDAC_Solver/HiDAC_PetscSolver.f90
> \
> # /Users/yaxiong/Documents/HiDAC_Code_Git/HiDAC_Tao_NURBS_PETSC/Hidac_Fortran_Lib/HiDACX/obj/debug/HiDAC_baseSolver.o
> \
> # /Users/yaxiong/Documents/HiDAC_Code_Git/HiDAC_Tao_NURBS_PETSC/Hidac_Fortran_Lib/HiDACX/obj/debug/HiDAC_Utils.o
> \
> #    This part is removed
>
> endif
>
> #============================Makefile
> finished ===============================
>
>
>
> and the main program starts like this:
>
>
> !============================  Main program starts=========================
>
> program main
>   !******************************************************!
>   ! Import all HiDAC library
>   !******************************************************!
>   include 'mpif.h'
> #include <petsc/finclude/petscksp.h>
>   use petscksp
>   use read_iges_file
>   use HiDAC_Approximation
>   use HiDAC_Domain
>   use HiDAC_Material
>   use HiDAC_Analysis
>   use HiDAC_Utils
>   use HiDAC_quadType
>
>   implicit none
> !====================*=============*================================
>
> I have not call the function in PETSC yet. Just add the KSP library in my
> main program. But I received  the following error when compiling:
>
>
> HiDACX/obj/debug/optimal_mechanical_part.o
> mpif90 -cpp -JHiDACX/obj/debug/  -Wall -fcheck=bounds -I/usr/local/include
> -I/Users/yaxiong/Downloads/petsc-3.9.4/include -c
> /Users/yaxiong/Documents/HiDAC_Code_Git/HiDAC_Tao_NURBS_Petsc/Hidac_Fortran_Lib/HiDAC_Test/optimal_mechanical_part.f90
> -o HiDACX/obj/debug/optimal_mechanical_part.o
>
> /Users/yaxiong/Downloads/petsc-3.9.4/include/petsc/finclude/petscsys.h:13:2:
>
>  #if defined (PETSC_HAVE_MPIUNI)
>   1~~~~~~~~~~~~
> Fatal Error: petscconf.h: No such file or directory
> compilation terminated.
> make[1]: *** [HiDACX/obj/debug/optimal_mechanical_part.o] Error 1
> make: *** [debug] Error 2
>
> How can I solve this problem?
>

You named the file

  optimal_mechanical_part.f90

but it has to be

  optimal_mechanical_part.F90

in order for the compiler to preprocess it.

   Matt


> Thanks
>
>
> Yaxiong Chen,
> Ph.D. Student
>
> School of Mechanical Engineering, 3171
>
> 585 Purdue Mall
>
> West Lafayette, IN 47907
>
>
>
>
>
> ------------------------------
> *From:* Matthew Knepley <knepley at gmail.com>
> *Sent:* Friday, January 18, 2019 2:21 PM
> *To:* Yaxiong Chen; PETSc
> *Subject:* Re: [petsc-users]Integrate PETSC with existing Fortran
>
> On Fri, Jan 18, 2019 at 2:06 PM Yaxiong Chen <chen2018 at purdue.edu> wrote:
>
> Hi Matt,
>
>
> Do not drop the Cc.
>
>
> I have hundreds .F90 files to compile as well as an .c file to compile in
> this project. The present makefile is more complex than the PETSC makefile.
> Do you mean I should use the PETSC makefile as a template and add all the
> dependencies in it?
>
>
> If you want to do it by hand you can, its just more error prone and
> verbose.
>
>
> My present Makefile is like the following(Red parts are news added related
> to Petsc). Is this the right direction? Since a static pattern is defined
> while the compiling for PetscSolver is a little different. I guess I should
> list it as a specific case instead of using the pattern.
>
>
> I will comment inline.
>
>
> #
> ============================================================================
> # Name        : Makefile
> # Author      : Tao Song
> # Version     : 2.1
> # Copyright   : GNU 2.1
> # Description : Makefile for HiDACX
> #
> ============================================================================
>
> OUT_DIR := HiDACX/
>
> ## Define mode
> ifndef MODE
>
> .PHONY: debug release
>
> debug: export MODE := debug
> release: export MODE := release
> debug release:
> @$(MAKE) $@
>
> clean:
> rm -rf $(OUT_DIR)obj/ $(OUT_DIR)lib/ $(OUT_DIR)bin/
>
> else
>
> ## SRCS = $(shell find ../../MyGit/Hidac_Fortran_Lib/ "*.f90")
> SRCS_DIR :=
> $(HOME)/Documents/HiDAC_Code_Git/HiDAC_Tao_NURBS_Petsc/Hidac_Fortran_Lib/
> FSRCS := \
> HiDAC_Utils/HiDAC_parameter.f90 \
> HiDAC_Utils/HiDAC_dataStructure.f90 \
> ......
> <other .f90 files>
> HiDAC_Domain/HiDAC_primitive.f90 \
> HiDAC_Solver/HiDAC_PetscSolver.F90 \
> HiDAC_Solver/HiDAC_Solver.f90 \
> IGES/read_iges_file.f90 \
> IGES/nrb.f90 \
> IGES/parse.f90 \
> IGES/map.f90 \
> DXF/dxf_reader.f90 \
>
> TSRCS_DIR := $(SRCS_DIR)HiDAC_Test/
>
> CSRCS := \
> HiDAC_Utils/HiDAC_DistanceField/HiDAC_QhullInterface_C.c \
>
> SRCS := $(addprefix $(SRCS_DIR), $(FSRCS) $(CSRCS))
>
> OBJS_DIR := $(OUT_DIR)obj/$(MODE)/
> LIB_DIR := $(OUT_DIR)lib/$(MODE)/
> EXE_DIR := $(OUT_DIR)bin/$(MODE)/
> EXT_DIR := /usr/local/
>
> #petsc lib
> include ${PETSC_DIR}/lib/petsc/conf/variables
> include ${PETSC_DIR}/lib/petsc/conf/rules
>
>
> You do not need 'rules' if you are doing things by hand. Although it is
> easier to just use the builtin rule for building source.
>
>
>
> #========================
>
> OBJS := $(addprefix $(OBJS_DIR), \
> $(patsubst %.f90, %.o, $(notdir $(FSRCS))) \
> $(patsubst %.c, %.o, $(notdir $(CSRCS))))
>
> LIB = $(LIB_DIR)libhidac.a
>
> VPATH := $(sort $(dir $(SRCS))) $(TSRCS_DIR)
>
> FC := mpif90
> CC := mpicc
> AR := ar rcs
>
> FCFLAGS := -cpp -J$(OBJS_DIR) # ${PETSC_KSP_LIB}# -fprofile-arcs
> -ftest-coverage
> CCFLAGS := -Dqh_QHpointer -ansi
> LFLAGS := -Wl,-rpath,$(EXT_DIR)lib, ${PETSC_KSP_LIB}
>
> FIDIR := \
> -I$(EXT_DIR)include \
> -I$(PETSC_DIR)/include
> CIDIR := \
> -I$(EXT_DIR)include/libqhull
>
> # sequence does matter
> EXTLIBS := \
> -L$(EXT_DIR)lib \
> -ldmumps -lmumps_common \
> -lparmetis -lmetis -lpord -lscalapack \
> -llapack -lblas -lqhullstatic_p \
> -lpetsc
>
>
> .PHONY: debug release all clean
>
> # Handle different mode using target-specific variable
> debug: FCFLAGS += -Wall -fcheck=bounds
> debug: CCFLAGS += -Wall
> debug: LFLAGS += -Wall -fcheck=bounds
>
> release: FCFLAGS += -O2
> release: CCFLAGS += -O2
> release: LFLAGS += -O2
>
> debug release: $(EXE_DIR)optimal_mechanical_part.hdc
>
> ## Rules
> .SUFFIXES: .hdc
> # Note the sequence of $(OBJS_DIR)%.o $(LIB) and $(EXTLIBS) does matter!
> Lib later!
> $(EXE_DIR)%.hdc: $(OBJS_DIR)%.o $(LIB)
> @mkdir -p $(EXE_DIR)
> echo $OBJS_DIR
> $(FC) $(LFLAGS) -o $@ $^ $(EXTLIBS)
>
> $(LIB): $(OBJS)
> echo $OBJS
> @mkdir -p $(LIB_DIR)
> $(AR) $@ $^
>
> $(OBJS_DIR)%.o:
> @mkdir -p $(OBJS_DIR)
> @echo $@
> $(if $(findstring .f90, $<), \
> $(FC) $(FCFLAGS) $(FIDIR) -c $< -o $@, \
> $(CC) $(CCFLAGS) $(CIDIR) -c $< -o $@)
>
> ## Dependencies of files
> $(OBJS_DIR)read_iges_file.o: \
> read_iges_file.f90 \
> $(OBJS_DIR)HiDAC_Approximation.o \
> $(OBJS_DIR)HiDAC_Utils.o \
> $(OBJS_DIR)nrb.o \
> $(OBJS_DIR)parse.o \
> $(OBJS_DIR)map.o \
> $(OBJS_DIR)dxf_reader.o
>
> $(OBJS_DIR)test_main.o: \
> test_main.f90 \
> $(OBJS_DIR)read_iges_file.o
>
> <.... Other dependencies>
>
> $ (OBJS_DIR)HiDAC_PetscSolver.o : HiDAC_PetscSolver.F90
> $(OBJS_DIR)HiDAC_baseSolver.o $(OBJS_DIR)HiDAC_Utils.o
> mpif90 -cpp -JHiDACX/obj/debug/ -Wall -ffree-line-length-0
> -Wno-unused-dummy-argument -g
> -I/Users/yaxiong/Downloads/petsc-3.9.4/include \
> -I/Users/yaxiong/Downloads/petsc-3.9.4/arch-darwin-c-debug/include \
> -I/opt/X11/include -o $(OBJS_DIR)HiDAC_PetscSolver.o \
> /Users/yaxiong/Documents/HiDAC_Code_Git/HiDAC_Tao_NURBS_PETSC/Hidac_Fortran_Lib/HiDAC_Solver/HiDAC_PetscSolver.F90
> \
> /Users/yaxiong/Documents/HiDAC_Code_Git/HiDAC_Tao_NURBS_PETSC/Hidac_Fortran_Lib/HiDACX/obj/debug/HiDAC_baseSolver.o
> \
> /Users/yaxiong/Documents/HiDAC_Code_Git/HiDAC_Tao_NURBS_PETSC/Hidac_Fortran_Lib/HiDACX/obj/debug/HiDAC_Utils.o
> \
> ${PETSC_KSP_LIB}
>
>
> You can do this, but to be certain you should should use
>
>   ${FC} -c ${FC_FLAGS} ${FFLAGS} ${FCPPFLAGS}
>
> along with your flags (-J, etc.). Notice that you do not need the *.o
> files here since you are only compiling.
>
> You can simplify this by includeing 'rules' and just putting the $
> (OBJS_DIR)HiDAC_PetscSolver.o as a
> dependence for your link.
>
>   Thanks,
>
>     Matt
>
>
> endif
>
>
> Yaxiong Chen,
> Ph.D. Student
>
> School of Mechanical Engineering, 3171
>
> 585 Purdue Mall
>
> West Lafayette, IN 47907
>
>
>
>
> ------------------------------
> *From:* Matthew Knepley <knepley at gmail.com>
> *Sent:* Friday, January 18, 2019 11:57 AM
> *To:* Yaxiong Chen; PETSc
> *Subject:* Re: [petsc-users]Integrate PETSC with existing Fortran
>
> On Fri, Jan 18, 2019 at 11:37 AM Yaxiong Chen <chen2018 at purdue.edu> wrote:
>
> Hi Matt,
>
> Now I have some problem with makefile.I want to integrate PETSC with a
> previous FEM Fortran code. To do that, I add a module named PetscSolver. It
> starts as the following:
>
>
> It needs to have suffix .F or .F90 so the preprocessor is invoked. Use the
> PETSc makefiles to avoid problems.
>
>   Matt
>
>
> !===================================================================
>
>   module HiDAC_PetscSolver
> #include <petsc/finclude/petscksp.h>
>   use petscksp
>   use HiDAC_Utils
>   use HiDAC_baseSolver
>     implicit none
>   private
>   type, extends(Solver), public  :: PetscSolver
>     private
>      Vec              xvec,bvec,uvec
>      Mat              Amat
>      KSP              ksp
>   contains
>      procedure :: initialize => initialize_PetscSolver
>      procedure :: assemble => assemble_PetscSolver
>      procedure :: solve => solve_PetscSolver
>      procedure :: postprocess => postprocess_PetscSolver
>      procedure :: finalize => finalize_PetscSolver
>      procedure :: getID
>   end type PetscSolver
> ......
> <Functions>
> end module HiDAC_PetscSolver
> !===================================================================
> However, it did not compile successfully. It gave me the error:
>
>  #include <petsc/finclude/petscksp.h>
>   1
> Warning: Illegal preprocessor directive
>
> /Users/yaxiong/Documents/HiDAC_Code_Git/HiDAC_Tao_NURBS_PETSC/Hidac_Fortran_Lib/HiDAC_Solver/HiDAC_PetscSolver.f90:23:5:
>
>       Vec              xvec,bvec,uvec
>      1
> Error: Unclassifiable statement at (1)
>
> /Users/yaxiong/Documents/HiDAC_Code_Git/HiDAC_Tao_NURBS_PETSC/Hidac_Fortran_Lib/HiDAC_Solver/HiDAC_PetscSolver.f90:24:5:
>
>       Mat              Amat
> All function and type related to Petsc are not defined.
>
> I must have made a very silly mistake. Can you show me how to fix it?
>
> Thanks
>
> Yaxiong Chen,
> Ph.D. Student
>
> School of Mechanical Engineering, 3171
>
> 585 Purdue Mall
>
> West Lafayette, IN 47907
>
>
>
>
> --
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> -- Norbert Wiener
>
> https://www.cse.buffalo.edu/~knepley/
> <http://www.cse.buffalo.edu/~knepley/>
>
>
>
> --
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> -- Norbert Wiener
>
> https://www.cse.buffalo.edu/~knepley/
> <http://www.cse.buffalo.edu/~knepley/>
>


-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener

https://www.cse.buffalo.edu/~knepley/ <http://www.cse.buffalo.edu/~knepley/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20190123/1de68e92/attachment-0001.html>


More information about the petsc-users mailing list